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