starMesh.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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "starMesh.H"
27 #include "emptyPolyPatch.H"
28 #include "demandDrivenData.H"
29 #include "cellModeller.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 // Merge tolerances
34 const Foam::scalar Foam::starMesh::smallMergeTol_ = 1e-3;
35 const Foam::scalar Foam::starMesh::cpMergePointTol_ = 1e-4;
36 
37 // Cell shape models
38 const Foam::cellModel* Foam::starMesh::unknownPtr_ =
39  Foam::cellModeller::lookup("unknown");
40 const Foam::cellModel* Foam::starMesh::tetPtr_ =
42 const Foam::cellModel* Foam::starMesh::pyrPtr_ =
44 const Foam::cellModel* Foam::starMesh::tetWedgePtr_ =
45  Foam::cellModeller::lookup("tetWedge");
46 const Foam::cellModel* Foam::starMesh::prismPtr_ =
48 const Foam::cellModel* Foam::starMesh::wedgePtr_ =
50 const Foam::cellModel* Foam::starMesh::hexPtr_ =
52 
53 const Foam::cellModel* Foam::starMesh::sammTrim1Ptr_ =
54  Foam::cellModeller::lookup("sammTrim1");
55 const Foam::cellModel* Foam::starMesh::sammTrim2Ptr_ =
56  Foam::cellModeller::lookup("sammTrim2");
57 const Foam::cellModel* Foam::starMesh::sammTrim3Ptr_ =
58  Foam::cellModeller::lookup("sammTrim3");
59 const Foam::cellModel* Foam::starMesh::sammTrim4Ptr_ =
60  Foam::cellModeller::lookup("sammTrim4");
61 const Foam::cellModel* Foam::starMesh::sammTrim5Ptr_ =
62  Foam::cellModeller::lookup("sammTrim5");
63 const Foam::cellModel* Foam::starMesh::sammTrim8Ptr_ =
64  Foam::cellModeller::lookup("hexagonalPrism");
65 
66 // Regular cell point addressing
67 // SAMM point addressing
68 const Foam::label Foam::starMesh::regularAddressingTable[6][8] =
69 {
70  { 0, 1, 2, 4, -1, -1, -1, -1}, // tet
71  { 0, 1, 2, 3, 4, -1, -1, -1}, // pyramid
72  { 0, 1, 2, 4, 6, -1, -1, -1}, // tet wedge
73  { 0, 1, 2, 4, 5, 6, -1, -1}, // prism
74  { 7, 6, 5, 3, 2, 1, 0, -1}, // wedge
75  { 0, 1, 2, 3, 4, 5, 6, 7} // hex
76 };
77 
78 
79 // SAMM point addressing
80 const Foam::label Foam::starMesh::sammAddressingTable[9][12] =
81 {
82  {-1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1}, // samm0 - empty
83  { 3, 2, 6, 7, 11, 9, 1, 5, 4, 12, -1, -1}, // samm1+
84  {13, 5, 6, 2, 10, 12, 4, 7, 3, 11, -1, -1}, // samm2+
85  { 2, 3, 0, 1, 10, 11, 12, 4, 8, 9, -1, -1}, // samm3+
86  { 0, 1, 3, 4, 13, 8, 9, 10, 11, 12, -1, -1}, // samm4+
87  {12, 7, 6, 5, 8, 11, 10, 9, -1, -1, -1, -1}, // samm5+
88  {-1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1}, // samm6 - empty
89  {-1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1}, // samm7 - empty
90  {11, 3, 15, 12, 4, 8, 10, 2, 14, 13, 5, 9} // samm8+
91 };
92 
93 
94 // lookup table giving OpenFOAM face number when looked up with shape index
95 // (first index) and STAR face number
96 // - first column is always -1
97 // - last column is -1 for all but hexagonal prism
98 // WARNING: Possible bug for sammTrim2
99 // The lookup table for SAMM shapes is based on the rotation of the
100 // shape. This would imply that the table below needs to be split between
101 // the regular shapes (3-9), which are OK, and the SAMM shapes, for which
102 // the face lookup needs to be done based on the rotation. Thus, for a samm
103 // cell, firts find out the face index in the normal rotation using the cell
104 // face permutation table and then use the index from the shape face lookup.
105 // Additionally, have in mind that this silliness does not allow matches
106 // on face 7 and 8 of the samm cell.
107 
108 const Foam::label Foam::starMesh::sammFacePermutationTable[24][8] =
109 {
110  {-1, 1, 2, 3, 4, 5, 6, 7}, // permutation 0
111  {-1, 3, 4, 5, 6, 1, 2, 7}, // permutation 1
112  {-1, 5, 6, 1, 2, 3, 4, 7}, // permutation 2
113  {-1, 1, 2, 5, 6, 4, 3, 7}, // permutation 3
114  {-1, 3, 4, 1, 2, 6, 5, 7}, // permutation 4
115  {-1, 5, 6, 3, 4, 2, 1, 7}, // permutation 5
116  {-1, 1, 2, 4, 3, 6, 5, 7}, // permutation 6
117  {-1, 3, 4, 6, 5, 2, 1, 7}, // permutation 7
118  {-1, 5, 6, 2, 1, 4, 3, 7}, // permutation 8
119  {-1, 1, 2, 6, 5, 3, 4, 7}, // permutation 9
120  {-1, 3, 4, 2, 1, 5, 6, 7}, // permutation 10
121  {-1, 5, 6, 4, 3, 1, 2, 7}, // permutation 11
122  {-1, 2, 1, 5, 6, 3, 4, 7}, // permutation 12
123  {-1, 4, 3, 1, 2, 5, 6, 7}, // permutation 13
124  {-1, 6, 5, 3, 4, 1, 2, 7}, // permutation 14
125  {-1, 2, 1, 3, 4, 6, 5, 7}, // permutation 15
126  {-1, 4, 3, 5, 6, 2, 1, 7}, // permutation 16
127  {-1, 6, 5, 1, 2, 4, 3, 7}, // permutation 17
128  {-1, 2, 1, 6, 5, 4, 3, 7}, // permutation 18
129  {-1, 4, 3, 2, 1, 6, 5, 7}, // permutation 19
130  {-1, 6, 5, 4, 3, 2, 1, 7}, // permutation 20
131  {-1, 2, 1, 4, 3, 5, 6, 7}, // permutation 21
132  {-1, 4, 3, 6, 5, 1, 2, 7}, // permutation 22
133  {-1, 6, 5, 2, 1, 3, 4, 7} // permutation 23
134 };
135 
136 const Foam::label Foam::starMesh::shapeFaceLookup[19][9] =
137 {
138  {-1, -1, -1, -1, -1, -1, -1, -1, -1}, // shape 0 - empty+
139  {-1, -1, -1, -1, -1, -1, -1, -1, -1}, // shape 1 - empty+
140  {-1, -1, -1, -1, -1, -1, -1, -1, -1}, // shape 2 - empty+
141  {-1, 4, 5, 2, 3, 0, 1, -1, -1}, // shape 3 - hex+
142  {-1, 1, 0, 5, 4, 2, 3, -1, -1}, // shape 4 - wedge+
143  {-1, 0, 1, 4, -1, 2, 3, -1, -1}, // shape 5 - prism+
144  {-1, 0, -1, 4, 2, 1, 3, -1, -1}, // shape 6 - pyr+
145  {-1, 3, -1, 2, -1, 1, 0, -1, -1}, // shape 7 - tet+
146  {-1, -1, -1, -1, -1, -1, -1, -1, -1}, // shape 8 - splitHex (empty)
147  {-1, 0, -1, 1, -1, 2, 3, -1, -1}, // shape 9 - tetWedge+
148  {-1, -1, -1, -1, -1, -1, -1, -1, -1}, // shape 10 - empty+
149  {-1, 1, 0, 3, 2, 5, 4, 6, -1}, // shape 11 - sammTrim1
150  {-1, 5, 4, 1, 0, 3, 2, 6, -1}, // shape 12 - sammTrim2
151  {-1, 2, 3, 0, 1, 4, 5, 6, -1}, // shape 13 - sammTrim3
152  {-1, 2, 3, 0, 1, 4, 5, 6, -1}, // shape 14 - sammTrim4
153  {-1, 5, 2, 4, 3, 1, 0, -1, -1}, // shape 15 - sammTrim5
154  {-1, -1, -1, -1, -1, -1, -1, -1, -1}, // shape 16 - empty
155  {-1, -1, -1, -1, -1, -1, -1, -1, -1}, // shape 17 - empty
156  {-1, 1, 0, 6, 7, 2, 3, 4, 5} // shape 18 - sammTrim8
157 };
158 
159 
160 // The star to foam face order mapping tables are potentially incomplete
161 // Currently available data is listed below.
162 // 1) hex and degenerate hex: OK
163 // samm trim 1:
164 // star number: 1 2 3 4 5 6 7 8 In ROTATION 0
165 // foam number: 5 4 1 0 3 2 6
166 // confirmed: 1 0 3 2 5 4 6
167 
168 // samm trim 2:
169 // star number: 1 2 3 4 5 6 7 8 In ROTATION 0
170 // foam number: 5 4 1 0 3 2 6
171 // confirmed: 4 0 3 2
172 
173 // samm trim 3:
174 // star number: 1 2 3 4 5 6 7 8 In ROTATION 0
175 // foam number: 2 3 0 1 4 5 6
176 // confirmed:
177 
178 // samm trim 4:
179 // star number: 1 2 3 4 5 6 7 8 In ROTATION 0
180 // foam number: 2 3 0 1 4 5 6
181 // confirmed:
182 
183 // samm trim 5:
184 // star number: 1 2 3 4 5 6 7 8 In ROTATION 0
185 // foam number: 2 4 3 1 0 5
186 // confirmed: 5 2 4 3 1 0
187 
188 // samm trim 8:
189 // star number: 1 2 3 4 5 6 7 8 In ROTATION 0
190 // foam number: 2 5 4 7 1 0 3 6
191 // confirmed: 1 0 6
192 
193 
194 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
195 
196 void Foam::starMesh::createPolyMeshData()
197 {
198  Info<< "Creating a polyMesh" << endl;
199 
200  createPolyCells();
201 
202  Info<< "\nNumber of internal faces: "
203  << nInternalFaces_ << endl;
204 
205  createPolyBoundary();
206 }
207 
208 
209 void Foam::starMesh::clearExtraStorage()
210 {
211  Info<< "Clearing extra storage" << endl;
212 
213  starPointLabelLookup_.setSize(0);
214  starPointID_.setSize(0);
215  starCellID_.setSize(0);
216  starCellLabelLookup_.setSize(0);
217  starCellPermutation_.setSize(0);
218  cellFaces_.setSize(0);
219  boundaryCellIDs_.setSize(0);
220  boundaryCellFaceIDs_.setSize(0);
221  couples_.clear();
222 
223  deleteDemandDrivenData(pointCellsPtr_);
224 }
225 
226 
227 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
228 
229 Foam::starMesh::starMesh
230 (
231  const fileName& prefix,
232  const Time& rt,
233  const scalar scaleFactor
234 )
235 :
236  casePrefix_(prefix),
237  runTime_(rt),
238  points_(0),
239  cellShapes_(0),
240  boundary_(0),
241  patchTypes_(0),
242  defaultFacesName_("defaultFaces"),
243  defaultFacesType_(emptyPolyPatch::typeName),
244  patchNames_(0),
245  patchPhysicalTypes_(0),
246  starPointLabelLookup_(0),
247  starPointID_(0),
248  starCellID_(0),
249  starCellLabelLookup_(0),
250  starCellPermutation_(0),
251  cellFaces_(0),
252  boundaryCellIDs_(0),
253  boundaryCellFaceIDs_(0),
254  meshFaces_(0),
255  cellPolys_(0),
256  nInternalFaces_(0),
257  polyBoundaryPatchStartIndices_(0),
258  pointCellsPtr_(NULL),
259  couples_(0),
260  isShapeMesh_(true)
261 {
262  readPoints(scaleFactor);
263 
264  readCells();
265 
266  readBoundary();
267 
268  fixCollapsedEdges();
269 
270  readCouples();
271 
272  if (couples_.size())
273  {
274  createCoupleMatches();
275  }
276 
277  markBoundaryFaces();
278 
279  mergeCoupleFacePoints();
280 
281  purgeCellShapes();
282 
283  collectBoundaryFaces();
284 }
285 
286 
287 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
288 
290 {
291  deleteDemandDrivenData(pointCellsPtr_);
292 }
293 
294 
295 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
static const cellModel * lookup(const word &)
Look up a model by name and return a pointer to the model or NULL.
Definition: cellModeller.C:88
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:78
~starMesh()
Destructor.
void setSize(const label)
Reset size of List.
Definition: List.C:295
Template functions to aid in the implementation of demand driven data.
Maps a geometry to a set of cell primitives, which enables geometric cell data to be calculated witho...
Definition: cellModel.H:64
messageStream Info
void deleteDemandDrivenData(DataPtr &dataPtr)