blockMeshCreate.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-2016 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 "error.H"
27 #include "blockMesh.H"
28 #include "cellModeller.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 void Foam::blockMesh::createPoints() const
33 {
34  const blockList& blocks = *this;
35 
36  if (verboseOutput)
37  {
38  Info<< "Creating points with scale " << scaleFactor_ << endl;
39  }
40 
41  //
42  // generate points
43  //
44  points_.clear();
45  points_.setSize(nPoints_);
46 
47  forAll(blocks, blockI)
48  {
49  const pointField& blockPoints = blocks[blockI].points();
50 
51  if (verboseOutput)
52  {
53  const Vector<label>& density = blocks[blockI].meshDensity();
54 
55  label v0 = blocks[blockI].vtxLabel(0, 0, 0);
56  label vi1 = blocks[blockI].vtxLabel(1, 0, 0);
57  scalar diStart = mag(blockPoints[vi1]-blockPoints[v0]);
58 
59  label vinM1 = blocks[blockI].vtxLabel(density.x()-1, 0, 0);
60  label vin = blocks[blockI].vtxLabel(density.x(), 0, 0);
61  scalar diFinal = mag(blockPoints[vin]-blockPoints[vinM1]);
62 
63  label vj1 = blocks[blockI].vtxLabel(0, 1, 0);
64  scalar djStart = mag(blockPoints[vj1]-blockPoints[v0]);
65  label vjnM1 = blocks[blockI].vtxLabel(0, density.y()-1, 0);
66  label vjn = blocks[blockI].vtxLabel(0, density.y(), 0);
67  scalar djFinal = mag(blockPoints[vjn]-blockPoints[vjnM1]);
68 
69  label vk1 = blocks[blockI].vtxLabel(0, 0, 1);
70  scalar dkStart = mag(blockPoints[vk1]-blockPoints[v0]);
71  label vknM1 = blocks[blockI].vtxLabel(0, 0, density.z()-1);
72  label vkn = blocks[blockI].vtxLabel(0, 0, density.z());
73  scalar dkFinal = mag(blockPoints[vkn]-blockPoints[vknM1]);
74 
75  Info<< " Block " << blockI << " cell size :" << nl
76  << " i : " << scaleFactor_*diStart << " .. "
77  << scaleFactor_*diFinal << nl
78  << " j : " << scaleFactor_*djStart << " .. "
79  << scaleFactor_*djFinal << nl
80  << " k : " << scaleFactor_*dkStart << " .. "
81  << scaleFactor_*dkFinal << nl
82  << endl;
83  }
84 
85  forAll(blockPoints, blockPointi)
86  {
87  points_
88  [
89  mergeList_
90  [
91  blockOffsets_[blockI] + blockPointi
92  ]
93  ] = scaleFactor_ * blockPoints[blockPointi];
94  }
95  }
96 }
97 
98 
99 void Foam::blockMesh::createCells() const
100 {
101  const blockList& blocks = *this;
102  const cellModel& hex = *(cellModeller::lookup("hex"));
103 
104  if (verboseOutput)
105  {
106  Info<< "Creating cells" << endl;
107  }
108 
109  //
110  // generate cells
111  //
112  cells_.clear();
113  cells_.setSize(nCells_);
114 
115  label cellLabel = 0;
116 
117  forAll(blocks, blockI)
118  {
119  const labelListList& blockCells = blocks[blockI].cells();
120 
121  forAll(blockCells, blockCelli)
122  {
123  labelList cellPoints(blockCells[blockCelli].size());
124 
125  forAll(cellPoints, cellPointi)
126  {
127  cellPoints[cellPointi] =
128  mergeList_
129  [
130  blockCells[blockCelli][cellPointi]
131  + blockOffsets_[blockI]
132  ];
133  }
134 
135  // Construct collapsed cell and add to list
136  cells_[cellLabel] = cellShape(hex, cellPoints, true);
137 
138  cellLabel++;
139  }
140  }
141 }
142 
143 
144 Foam::faceList Foam::blockMesh::createPatchFaces
145 (
146  const polyPatch& patchTopologyFaces
147 ) const
148 {
149  const blockList& blocks = *this;
150 
151  labelList blockLabels = patchTopologyFaces.polyPatch::faceCells();
152 
153  label nFaces = 0;
154 
155  forAll(patchTopologyFaces, patchTopologyFaceLabel)
156  {
157  const label blockI = blockLabels[patchTopologyFaceLabel];
158 
159  faceList blockFaces = blocks[blockI].blockShape().faces();
160 
161  forAll(blockFaces, blockFaceLabel)
162  {
163  if
164  (
165  blockFaces[blockFaceLabel]
166  == patchTopologyFaces[patchTopologyFaceLabel]
167  )
168  {
169  nFaces +=
170  blocks[blockI].boundaryPatches()[blockFaceLabel].size();
171  }
172  }
173  }
174 
175 
176  faceList patchFaces(nFaces);
177  face quadFace(4);
178  label faceLabel = 0;
179 
180  forAll(patchTopologyFaces, patchTopologyFaceLabel)
181  {
182  const label blockI = blockLabels[patchTopologyFaceLabel];
183 
184  faceList blockFaces = blocks[blockI].blockShape().faces();
185 
186  forAll(blockFaces, blockFaceLabel)
187  {
188  if
189  (
190  blockFaces[blockFaceLabel]
191  == patchTopologyFaces[patchTopologyFaceLabel]
192  )
193  {
194  const labelListList& blockPatchFaces =
195  blocks[blockI].boundaryPatches()[blockFaceLabel];
196 
197  forAll(blockPatchFaces, blockFaceLabel)
198  {
199  // Lookup the face points
200  // and collapse duplicate point labels
201 
202  quadFace[0] =
203  mergeList_
204  [
205  blockPatchFaces[blockFaceLabel][0]
206  + blockOffsets_[blockI]
207  ];
208 
209  label nUnique = 1;
210 
211  for
212  (
213  label facePointLabel = 1;
214  facePointLabel < 4;
215  facePointLabel++
216  )
217  {
218  quadFace[nUnique] =
219  mergeList_
220  [
221  blockPatchFaces[blockFaceLabel][facePointLabel]
222  + blockOffsets_[blockI]
223  ];
224 
225  if (quadFace[nUnique] != quadFace[nUnique-1])
226  {
227  nUnique++;
228  }
229  }
230 
231  if (quadFace[nUnique-1] == quadFace[0])
232  {
233  nUnique--;
234  }
235 
236  if (nUnique == 4)
237  {
238  patchFaces[faceLabel++] = quadFace;
239  }
240  else if (nUnique == 3)
241  {
242  patchFaces[faceLabel++] = face
243  (
245  );
246  }
247  // else the face has collapsed to an edge or point
248  }
249  }
250  }
251  }
252 
253  patchFaces.setSize(faceLabel);
254 
255  return patchFaces;
256 }
257 
258 
259 void Foam::blockMesh::createPatches() const
260 {
261  const polyPatchList& topoPatches = topology().boundaryMesh();
262 
263  if (verboseOutput)
264  {
265  Info<< "Creating patches" << endl;
266  }
267 
268  //
269  // generate points
270  //
271 
272  patches_.clear();
273  patches_.setSize(topoPatches.size());
274 
275  forAll(topoPatches, patchi)
276  {
277  patches_[patchi] = createPatchFaces(topoPatches[patchi]);
278  }
279 
280 }
281 
282 
284 {
285  blockList& blocks = *this;
286 
287  forAll(blocks, blockI)
288  {
289  blocks[blockI].clearGeom();
290  }
291 }
292 
293 // ************************************************************************* //
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:45
const polyMesh & topology() const
Return the blockMesh topology as a polyMesh.
Definition: blockMesh.C:82
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
List< face > faceList
Definition: faceListFwd.H:43
PtrList< block > blockList
A PtrList of blocks.
Definition: blockList.H:45
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
SubList< label > subList
Declare type of subList.
Definition: List.H:155
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:356
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const char nl
Definition: Ostream.H:262
face quadFace(4)
void setSize(const label)
Reset size of List.
Definition: List.C:295
label patchi
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
void clearGeom()
Clear geometry (internal points, cells, boundaryPatches)
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29