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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "blockMesh.H"
30 #include "cellModel.H"
31 #include "emptyPolyPatch.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void Foam::blockMesh::createPoints() const
36 {
37  const blockList& blocks = *this;
38 
39  if (verbose_)
40  {
41  Info<< "Creating points with scale " << scaleFactor_ << endl;
42  }
43 
44  points_.setSize(nPoints_);
45 
46  forAll(blocks, blocki)
47  {
48  const pointField& blockPoints = blocks[blocki].points();
49 
50  if (verbose_)
51  {
52  const label nx = blocks[blocki].density().x();
53  const label ny = blocks[blocki].density().y();
54  const label nz = blocks[blocki].density().z();
55 
56  label v0 = blocks[blocki].pointLabel(0, 0, 0);
57  label vi1 = blocks[blocki].pointLabel(1, 0, 0);
58  scalar diStart = mag(blockPoints[vi1] - blockPoints[v0]);
59 
60  label vinM1 = blocks[blocki].pointLabel(nx-1, 0, 0);
61  label vin = blocks[blocki].pointLabel(nx, 0, 0);
62  scalar diFinal = mag(blockPoints[vin] - blockPoints[vinM1]);
63 
64  label vj1 = blocks[blocki].pointLabel(0, 1, 0);
65  scalar djStart = mag(blockPoints[vj1] - blockPoints[v0]);
66  label vjnM1 = blocks[blocki].pointLabel(0, ny-1, 0);
67  label vjn = blocks[blocki].pointLabel(0, ny, 0);
68  scalar djFinal = mag(blockPoints[vjn] - blockPoints[vjnM1]);
69 
70  label vk1 = blocks[blocki].pointLabel(0, 0, 1);
71  scalar dkStart = mag(blockPoints[vk1] - blockPoints[v0]);
72  label vknM1 = blocks[blocki].pointLabel(0, 0, nz-1);
73  label vkn = blocks[blocki].pointLabel(0, 0, nz);
74  scalar dkFinal = mag(blockPoints[vkn] - blockPoints[vknM1]);
75 
76  Info<< " Block " << blocki << " cell size :" << nl
77  << " i : " << scaleFactor_*diStart << " .. "
78  << scaleFactor_*diFinal << nl
79  << " j : " << scaleFactor_*djStart << " .. "
80  << scaleFactor_*djFinal << nl
81  << " k : " << scaleFactor_*dkStart << " .. "
82  << scaleFactor_*dkFinal << nl
83  << endl;
84  }
85 
86  forAll(blockPoints, blockPointi)
87  {
88  points_
89  [
90  mergeList_
91  [
92  blockOffsets_[blocki] + blockPointi
93  ]
94  ] = scaleFactor_ * blockPoints[blockPointi];
95  }
96  }
97 }
98 
99 
100 void Foam::blockMesh::createCells() const
101 {
102  const blockList& blocks = *this;
103  const cellModel& hex = cellModel::ref(cellModel::HEX);
104 
105  if (verbose_)
106  {
107  Info<< "Creating cells" << endl;
108  }
109 
110  cells_.setSize(nCells_);
111 
112  label celli = 0;
113 
114  forAll(blocks, blocki)
115  {
116  const List<FixedList<label, 8>>& blockCells = blocks[blocki].cells();
117 
118  forAll(blockCells, blockCelli)
119  {
120  labelList cellPoints(blockCells[blockCelli].size());
121 
122  forAll(cellPoints, cellPointi)
123  {
124  cellPoints[cellPointi] =
125  mergeList_
126  [
127  blockCells[blockCelli][cellPointi]
128  + blockOffsets_[blocki]
129  ];
130  }
131 
132  // Construct collapsed cell and add to list
133  cells_[celli] = cellShape(hex, cellPoints, true);
134 
135  ++celli;
136  }
137  }
138 }
139 
140 
141 Foam::faceList Foam::blockMesh::createPatchFaces
142 (
143  const polyPatch& patchTopologyFaces
144 ) const
145 {
146  const blockList& blocks = *this;
147 
148  labelList blockLabels = patchTopologyFaces.polyPatch::faceCells();
149 
150  label nFaces = 0;
151 
152  forAll(patchTopologyFaces, patchTopologyFaceLabel)
153  {
154  const label blocki = blockLabels[patchTopologyFaceLabel];
155 
156  faceList blockFaces = blocks[blocki].blockShape().faces();
157 
158  forAll(blockFaces, blockFaceLabel)
159  {
160  if
161  (
162  blockFaces[blockFaceLabel]
163  == patchTopologyFaces[patchTopologyFaceLabel]
164  )
165  {
166  nFaces +=
167  blocks[blocki].boundaryPatches()[blockFaceLabel].size();
168  }
169  }
170  }
171 
172 
173  faceList patchFaces(nFaces);
174  face quadFace(4);
175  label faceLabel = 0;
176 
177  forAll(patchTopologyFaces, patchTopologyFaceLabel)
178  {
179  const label blocki = blockLabels[patchTopologyFaceLabel];
180 
181  faceList blockFaces = blocks[blocki].blockShape().faces();
182 
183  forAll(blockFaces, blockFaceLabel)
184  {
185  if
186  (
187  blockFaces[blockFaceLabel]
188  == patchTopologyFaces[patchTopologyFaceLabel]
189  )
190  {
191  const List<FixedList<label, 4>>& blockPatchFaces =
192  blocks[blocki].boundaryPatches()[blockFaceLabel];
193 
194  forAll(blockPatchFaces, blockFaceLabel)
195  {
196  // Lookup the face points
197  // and collapse duplicate point labels
198 
199  quadFace[0] =
200  mergeList_
201  [
202  blockPatchFaces[blockFaceLabel][0]
203  + blockOffsets_[blocki]
204  ];
205 
206  label nUnique = 1;
207 
208  for
209  (
210  label facePointLabel = 1;
211  facePointLabel < 4;
212  facePointLabel++
213  )
214  {
215  quadFace[nUnique] =
216  mergeList_
217  [
218  blockPatchFaces[blockFaceLabel][facePointLabel]
219  + blockOffsets_[blocki]
220  ];
221 
222  if (quadFace[nUnique] != quadFace[nUnique-1])
223  {
224  nUnique++;
225  }
226  }
227 
228  if (quadFace[nUnique-1] == quadFace[0])
229  {
230  nUnique--;
231  }
232 
233  if (nUnique == 4)
234  {
235  patchFaces[faceLabel++] = quadFace;
236  }
237  else if (nUnique == 3)
238  {
239  patchFaces[faceLabel++] = face
240  (
242  );
243  }
244  // else the face has collapsed to an edge or point
245  }
246  }
247  }
248  }
249 
250  patchFaces.setSize(faceLabel);
251 
252  return patchFaces;
253 }
254 
255 
256 void Foam::blockMesh::createPatches() const
257 {
258  const polyPatchList& topoPatches = topology().boundaryMesh();
259 
260  if (verbose_)
261  {
262  Info<< "Creating patches" << endl;
263  }
264 
265  patches_.setSize(topoPatches.size());
266 
267  forAll(topoPatches, patchi)
268  {
269  patches_[patchi] = createPatchFaces(topoPatches[patchi]);
270  }
271 }
272 
273 
274 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
275 
278 {
279  const blockMesh& blkMesh = *this;
280 
281  if (verbose_)
282  {
283  Info<< nl << "Creating polyMesh from blockMesh" << endl;
284  }
285 
287  (
288  io,
289  pointField(blkMesh.points()), // Copy, could we re-use space?
290  blkMesh.cells(),
291  blkMesh.patches(),
292  blkMesh.patchNames(),
293  blkMesh.patchDicts(),
294  "defaultFaces", // Default patch name
295  emptyPolyPatch::typeName // Default patch type
296  );
297 
298 
299  // Set any cellZones
300  const label nZones = blkMesh.numZonedBlocks();
301 
302  if (nZones)
303  {
304  polyMesh& pmesh = *meshPtr;
305 
306  if (verbose_)
307  {
308  Info<< "Adding cell zones" << endl;
309  }
310 
311  // Map from zoneName to cellZone index
312  HashTable<label> zoneMap(2*nZones);
313 
314  // Cells per zone
315  List<DynamicList<label>> zoneCells(nZones);
316 
317  // Running cell counter
318  label celli = 0;
319 
320  // Largest zone so far
321  label freeZonei = 0;
322 
323  for (const block& b : blkMesh)
324  {
325  const word& zoneName = b.zoneName();
326  const label nCellsInBlock = b.cells().size();
327 
328  if (zoneName.size())
329  {
330  const auto iter = zoneMap.cfind(zoneName);
331 
332  label zonei = freeZonei;
333 
334  if (iter.found())
335  {
336  zonei = *iter;
337  }
338  else
339  {
340  zoneMap.insert(zoneName, zonei);
341  ++freeZonei;
342 
343  if (verbose_)
344  {
345  Info<< " " << zonei << '\t' << zoneName << endl;
346  }
347  }
348 
349 
350  // Fill with cell ids
351 
352  zoneCells[zonei].reserve
353  (
354  zoneCells[zonei].size() + nCellsInBlock
355  );
356 
357  const label endOfFill = celli + nCellsInBlock;
358 
359  for (; celli < endOfFill; ++celli)
360  {
361  zoneCells[zonei].append(celli);
362  }
363  }
364  else
365  {
366  celli += nCellsInBlock;
367  }
368  }
369 
370  List<cellZone*> cz(zoneMap.size());
371  forAllConstIters(zoneMap, iter)
372  {
373  const word& zoneName = iter.key();
374  const label zonei = iter.val();
375 
376  cz[zonei] = new cellZone
377  (
378  zoneName,
379  zoneCells[zonei].shrink(),
380  zonei,
381  pmesh.cellZones()
382  );
383  }
384 
385  pmesh.pointZones().resize(0);
386  pmesh.faceZones().resize(0);
387  pmesh.cellZones().resize(0);
388  pmesh.addZones(List<pointZone*>(), List<faceZone*>(), cz);
389  }
390 
391 
392  // Merge patch pairs, cyclic must be done elsewhere
393  // - requires libdynamicMesh
394 
395  return meshPtr;
396 }
397 
398 
399 // ************************************************************************* //
nZones
label nZones
Definition: interpolatedFaces.H:24
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::HashTable::size
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Foam::blockMesh::patchDicts
PtrList< dictionary > patchDicts() const
Get patch information from the topology mesh.
Definition: blockMesh.C:153
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::block
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:58
meshPtr
Foam::autoPtr< Foam::fvMesh > meshPtr(nullptr)
Foam::blockMesh::mesh
autoPtr< polyMesh > mesh(const IOobject &io) const
Create polyMesh, with cell zones.
Definition: blockMeshCreate.C:277
Foam::blockMesh::patches
const faceListList & patches() const
Return the patch face lists.
Definition: blockMesh.C:198
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::cellModel::HEX
hex
Definition: cellModel.H:81
quadFace
face quadFace(4)
Foam::HashTable::insert
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:492
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
blockMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:62
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:486
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
cellModel.H
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyMesh::pointZones
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:480
Foam::polyPatchList
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:47
Foam::cellModel::ref
static const cellModel & ref(const modelType model)
Look up reference to cellModel by enumeration. Fatal on failure.
Definition: cellModels.C:157
Foam::List< label >::subList
SubList< label > subList
Declare type of subList.
Definition: List.H:114
emptyPolyPatch.H
Foam::hex
IOstream & hex(IOstream &io)
Definition: IOstream.H:437
Foam::HashTable< label >
Foam::blockMesh::blockList
PtrList< block > blockList
The list of blocks is stored as a PtrList.
Definition: blockMesh.H:157
Foam::autoPtr< Foam::polyMesh >
Foam::blockMesh::numZonedBlocks
label numZonedBlocks() const
Number of blocks with specified zones.
Definition: blockMesh.C:227
Foam::HashTable::cfind
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition: HashTableI.H:141
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::blockMesh::points
const pointField & points() const
The points for the entire mesh.
Definition: blockMesh.C:176
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::blockMesh::patchNames
wordList patchNames() const
Return patch names.
Definition: blockMesh.C:209
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::blockList
PtrList< block > blockList
A PtrList of blocks.
Definition: blockList.H:47
Foam::blockMesh::cells
const cellShapeList & cells() const
Return cell shapes list.
Definition: blockMesh.C:187
Foam::polyMesh::addZones
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:999
Foam::blockMesh
A multi-block mesh generator.
Definition: blockMesh.H:148