Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / applications / utilities / mesh / conversion / star3ToFoam / starMesh.C
blobbe5a147166e4a4feb9f1fcee4ec11b8d5b42a706
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 OpenCFD Ltd.
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 "starMesh.H"
27 #include "emptyPolyPatch.H"
28 #include "demandDrivenData.H"
29 #include "cellModeller.H"
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 // Merge tolerances
34 // Moni, tolerances:
35 const scalar starMesh::smallMergeTol_ = 1e-3;
36 const scalar starMesh::cpMergePointTol_ = 1e-4;
38 // const scalar starMesh::smallMergeTol_ = 4e-4;
39 // const scalar starMesh::cpMergePointTol_ = 1e-3;
41 // Cell shape models
42 const cellModel* starMesh::unknownPtr_ = cellModeller::lookup("unknown");
43 const cellModel* starMesh::tetPtr_ = cellModeller::lookup("tet");
44 const cellModel* starMesh::pyrPtr_ = cellModeller::lookup("pyr");
45 const cellModel* starMesh::tetWedgePtr_ = cellModeller::lookup("tetWedge");
46 const cellModel* starMesh::prismPtr_ = cellModeller::lookup("prism");
47 const cellModel* starMesh::wedgePtr_ = cellModeller::lookup("wedge");
48 const cellModel* starMesh::hexPtr_ = cellModeller::lookup("hex");
50 const cellModel* starMesh::sammTrim1Ptr_ = cellModeller::lookup("sammTrim1");
51 const cellModel* starMesh::sammTrim2Ptr_ = cellModeller::lookup("sammTrim2");
52 const cellModel* starMesh::sammTrim3Ptr_ = cellModeller::lookup("sammTrim3");
53 const cellModel* starMesh::sammTrim4Ptr_ = cellModeller::lookup("sammTrim4");
54 const cellModel* starMesh::sammTrim5Ptr_ = cellModeller::lookup("sammTrim5");
55 const cellModel* starMesh::sammTrim8Ptr_ =
56     cellModeller::lookup("hexagonalPrism");
58 // Regular cell point addressing
59 // SAMM point addressing
60 const label starMesh::regularAddressingTable[6][8] =
62     { 0,  1,  2,  4, -1, -1, -1, -1},    // tet
63     { 0,  1,  2,  3,  4, -1, -1, -1},    // pyramid
64     { 0,  1,  2,  4,  6, -1, -1, -1},    // tet wedge
65     { 0,  1,  2,  4,  5,  6, -1, -1},    // prism
66     { 7,  6,  5,  3,  2,  1,  0, -1},    // wedge
67     { 0,  1,  2,  3,  4,  5,  6,  7}     // hex
71 // SAMM point addressing
72 const label starMesh::sammAddressingTable[9][12] =
74     {-1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1},    // samm0 - empty
75     { 3,  2,  6,  7, 11,  9,  1,  5,  4, 12, -1, -1},    // samm1+
76     {13,  5,  6,  2, 10, 12,  4,  7,  3, 11, -1, -1},    // samm2+
77     { 2,  3,  0,  1, 10, 11, 12,  4,  8,  9, -1, -1},    // samm3+
78     { 0,  1,  3,  4, 13,  8,  9, 10, 11, 12, -1, -1},    // samm4+
79     {12,  7,  6,  5,  8, 11, 10,  9, -1, -1, -1, -1},    // samm5+
80     {-1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1},    // samm6 - empty
81     {-1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1},    // samm7 - empty
82     {11,  3, 15, 12,  4,  8, 10,  2, 14, 13,  5,  9}     // samm8+
86 // lookup table giving OpenFOAM face number when looked up with shape index
87 // (first index) and STAR face number
88 // - first column is always -1
89 // - last column is -1 for all but hexagonal prism
90 // WARNING: Possible bug for sammTrim2
91 // The lookup table for SAMM shapes is based on the rotation of the
92 // shape. This would imply that the table below needs to be split between
93 // the regular shapes (3-9), which are OK, and the SAMM shapes, for which
94 // the face lookup needs to be done based on the rotation. Thus, for a samm
95 // cell, firts find out the face index in the normal rotation using the cell
96 // face permutation table and then use the index from the shape face lookup.
97 // Additionally, have in mind that this silliness does not allow matches
98 // on face 7 and 8 of the samm cell.
100 const label starMesh::sammFacePermutationTable[24][8] =
102   {-1, 1, 2, 3, 4, 5, 6, 7},    // permutation   0
103   {-1, 3, 4, 5, 6, 1, 2, 7},    // permutation   1
104   {-1, 5, 6, 1, 2, 3, 4, 7},    // permutation   2
105   {-1, 1, 2, 5, 6, 4, 3, 7},    // permutation   3
106   {-1, 3, 4, 1, 2, 6, 5, 7},    // permutation   4
107   {-1, 5, 6, 3, 4, 2, 1, 7},    // permutation   5
108   {-1, 1, 2, 4, 3, 6, 5, 7},    // permutation   6
109   {-1, 3, 4, 6, 5, 2, 1, 7},    // permutation   7
110   {-1, 5, 6, 2, 1, 4, 3, 7},    // permutation   8
111   {-1, 1, 2, 6, 5, 3, 4, 7},    // permutation   9
112   {-1, 3, 4, 2, 1, 5, 6, 7},    // permutation  10
113   {-1, 5, 6, 4, 3, 1, 2, 7},    // permutation  11
114   {-1, 2, 1, 5, 6, 3, 4, 7},    // permutation  12
115   {-1, 4, 3, 1, 2, 5, 6, 7},    // permutation  13
116   {-1, 6, 5, 3, 4, 1, 2, 7},    // permutation  14
117   {-1, 2, 1, 3, 4, 6, 5, 7},    // permutation  15
118   {-1, 4, 3, 5, 6, 2, 1, 7},    // permutation  16
119   {-1, 6, 5, 1, 2, 4, 3, 7},    // permutation  17
120   {-1, 2, 1, 6, 5, 4, 3, 7},    // permutation  18
121   {-1, 4, 3, 2, 1, 6, 5, 7},    // permutation  19
122   {-1, 6, 5, 4, 3, 2, 1, 7},    // permutation  20
123   {-1, 2, 1, 4, 3, 5, 6, 7},    // permutation  21
124   {-1, 4, 3, 6, 5, 1, 2, 7},    // permutation  22
125   {-1, 6, 5, 2, 1, 3, 4, 7}     // permutation  23
128 const label starMesh::shapeFaceLookup[19][9] =
130     {-1, -1, -1, -1, -1, -1, -1, -1, -1},    // shape  0 - empty+
131     {-1, -1, -1, -1, -1, -1, -1, -1, -1},    // shape  1 - empty+
132     {-1, -1, -1, -1, -1, -1, -1, -1, -1},    // shape  2 - empty+
133     {-1,  4,  5,  2,  3,  0,  1, -1, -1},    // shape  3 - hex+
134     {-1,  1,  0,  5,  4,  2,  3, -1, -1},    // shape  4 - wedge+
135     {-1,  0,  1,  4, -1,  2,  3, -1, -1},    // shape  5 - prism+
136     {-1,  0, -1,  4,  2,  1,  3, -1, -1},    // shape  6 - pyr+
137     {-1,  3, -1,  2, -1,  1,  0, -1, -1},    // shape  7 - tet+
138     {-1, -1, -1, -1, -1, -1, -1, -1, -1},    // shape  8 - splitHex (empty)
139     {-1,  0, -1,  1, -1,  2,  3, -1, -1},    // shape  9 - tetWedge+
140     {-1, -1, -1, -1, -1, -1, -1, -1, -1},    // shape 10 - empty+
141     {-1,  1,  0,  3,  2,  5,  4,  6, -1},    // shape 11 - sammTrim1
142     {-1,  5,  4,  1,  0,  3,  2,  6, -1},    // shape 12 - sammTrim2
143     {-1,  2,  3,  0,  1,  4,  5,  6, -1},    // shape 13 - sammTrim3
144     {-1,  2,  3,  0,  1,  4,  5,  6, -1},    // shape 14 - sammTrim4
145     {-1,  5,  2,  4,  3,  1,  0, -1, -1},    // shape 15 - sammTrim5
146     {-1, -1, -1, -1, -1, -1, -1, -1, -1},    // shape 16 - empty
147     {-1, -1, -1, -1, -1, -1, -1, -1, -1},    // shape 17 - empty
148     {-1,  1,  0,  6,  7,  2,  3,  4,  5}     // shape 18 - sammTrim8
152 // The star to foam face order mapping tables are potentially incomplete
153 // Currently available data is listed below.
154 // 1) hex and degenerate hex: OK
155 // samm trim 1:
156 // star number: 1 2 3 4 5 6 7 8  In ROTATION 0
157 // foam number: 5 4 1 0 3 2 6
158 // confirmed:   1 0 3 2 5 4 6
160 // samm trim 2:
161 // star number: 1 2 3 4 5 6 7 8  In ROTATION 0
162 // foam number: 5 4 1 0 3 2 6
163 // confirmed:     4   0 3 2
165 // samm trim 3:
166 // star number: 1 2 3 4 5 6 7 8  In ROTATION 0
167 // foam number: 2 3 0 1 4 5 6
168 // confirmed:
170 // samm trim 4:
171 // star number: 1 2 3 4 5 6 7 8  In ROTATION 0
172 // foam number: 2 3 0 1 4 5 6
173 // confirmed:
175 // samm trim 5:
176 // star number: 1 2 3 4 5 6 7 8  In ROTATION 0
177 // foam number: 2 4 3 1 0 5
178 // confirmed:   5 2 4 3 1 0
180 // samm trim 8:
181 // star number: 1 2 3 4 5 6 7 8  In ROTATION 0
182 // foam number: 2 5 4 7 1 0 3 6
183 // confirmed:   1 0 6
186 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
188 // Make polyhedral mesh data (packing)
189 void starMesh::createPolyMeshData()
191     Info<< "Creating a polyMesh" << endl;
193     createPolyCells();
195     Info<< "\nNumber of internal faces: "
196         << nInternalFaces_ << endl;
198     createPolyBoundary();
202 // Clear extra storage before creation of the mesh to remove
203 // a memory peak
204 void starMesh::clearExtraStorage()
206     Info<< "Clearing extra storage" << endl;
208     starPointLabelLookup_.setSize(0);
209     starPointID_.setSize(0);
210     starCellID_.setSize(0);
211     starCellLabelLookup_.setSize(0);
212     starCellPermutation_.setSize(0);
213     cellFaces_.setSize(0);
214     boundaryCellIDs_.setSize(0);
215     boundaryCellFaceIDs_.setSize(0);
216     couples_.clear();
218     deleteDemandDrivenData(pointCellsPtr_);
221 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
223 // Construct from components
224 starMesh::starMesh
226     const fileName& prefix,
227     const Time& rt,
228     const scalar scaleFactor
231     casePrefix_(prefix),
232     runTime_(rt),
233     points_(0),
234     cellShapes_(0),
235     boundary_(0),
236     patchTypes_(0),
237     defaultFacesName_("defaultFaces"),
238     defaultFacesType_(emptyPolyPatch::typeName),
239     patchNames_(0),
240     patchPhysicalTypes_(0),
241     starPointLabelLookup_(0),
242     starPointID_(0),
243     starCellID_(0),
244     starCellLabelLookup_(0),
245     starCellPermutation_(0),
246     cellFaces_(0),
247     boundaryCellIDs_(0),
248     boundaryCellFaceIDs_(0),
249     meshFaces_(0),
250     cellPolys_(0),
251     nInternalFaces_(0),
252     polyBoundaryPatchStartIndices_(0),
253     pointCellsPtr_(NULL),
254     couples_(0),
255     isShapeMesh_(true)
257     readPoints(scaleFactor);
259     readCells();
261     readBoundary();
263     fixCollapsedEdges();
265     readCouples();
267     if (couples_.size())
268     {
269         createCoupleMatches();
270     }
272     markBoundaryFaces();
274     mergeCoupleFacePoints();
276     purgeCellShapes();
278     collectBoundaryFaces();
282 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
284 starMesh::~starMesh()
286     deleteDemandDrivenData(pointCellsPtr_);
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 // ************************************************************************* //