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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
35void Foam::blockMesh::createPoints() const
36{
37 const blockList& blocks = *this;
38
39 const vector scaleTot
40 (
41 prescaling_.x() * scaling_.x(),
42 prescaling_.y() * scaling_.y(),
43 prescaling_.z() * scaling_.z()
44 );
45
46 if (verbose_)
47 {
48 Info<< "Creating points with scale " << scaleTot << endl;
49 }
50
51 points_.resize(nPoints_);
52
53 forAll(blocks, blocki)
54 {
55 const pointField& blockPoints = blocks[blocki].points();
56
57 const labelSubList pointAddr
58 (
59 mergeList_,
60 blockPoints.size(),
61 blockOffsets_[blocki]
62 );
63
64 UIndirectList<point>(points_, pointAddr) = blockPoints;
65
66 if (verbose_)
67 {
68 Info<< " Block " << blocki << " cell size :" << nl;
69
70 const label v0 = blocks[blocki].pointLabel(0, 0, 0);
71
72 {
73 const label nx = blocks[blocki].density().x();
74 const label v1 = blocks[blocki].pointLabel(1, 0, 0);
75 const label vn = blocks[blocki].pointLabel(nx, 0, 0);
76 const label vn1 = blocks[blocki].pointLabel(nx-1, 0, 0);
77
78 const scalar cwBeg = mag(blockPoints[v1] - blockPoints[v0]);
79 const scalar cwEnd = mag(blockPoints[vn] - blockPoints[vn1]);
80
81 Info<< " i : "
82 << cwBeg*scaleTot.x() << " .. "
83 << cwEnd*scaleTot.x() << nl;
84 }
85
86 {
87 const label ny = blocks[blocki].density().y();
88 const label v1 = blocks[blocki].pointLabel(0, 1, 0);
89 const label vn = blocks[blocki].pointLabel(0, ny, 0);
90 const label vn1 = blocks[blocki].pointLabel(0, ny-1, 0);
91
92 const scalar cwBeg = mag(blockPoints[v1] - blockPoints[v0]);
93 const scalar cwEnd = mag(blockPoints[vn] - blockPoints[vn1]);
94
95 Info<< " j : "
96 << cwBeg*scaleTot.y() << " .. "
97 << cwEnd*scaleTot.y() << nl;
98 }
99
100 {
101 const label nz = blocks[blocki].density().z();
102 const label v1 = blocks[blocki].pointLabel(0, 0, 1);
103 const label vn = blocks[blocki].pointLabel(0, 0, nz);
104 const label vn1 = blocks[blocki].pointLabel(0, 0, nz-1);
105
106 const scalar cwBeg = mag(blockPoints[v1] - blockPoints[v0]);
107 const scalar cwEnd = mag(blockPoints[vn] - blockPoints[vn1]);
108
109 Info<< " k : "
110 << cwBeg*scaleTot.z() << " .. "
111 << cwEnd*scaleTot.z() << nl;
112 }
113 Info<< endl;
114 }
115 }
116
117 inplacePointTransforms(points_);
118}
119
120
121void Foam::blockMesh::createCells() const
122{
123 const blockList& blocks = *this;
124 const cellModel& hex = cellModel::ref(cellModel::HEX);
125
126 if (verbose_)
127 {
128 Info<< "Creating cells" << endl;
129 }
130
131 cells_.resize(nCells_);
132
133 label celli = 0;
134
135 labelList cellPoints(8); // Hex cells - 8 points
136
137 forAll(blocks, blocki)
138 {
139 for (const hexCell& blockCell : blocks[blocki].cells())
140 {
141 forAll(cellPoints, cellPointi)
142 {
143 cellPoints[cellPointi] =
144 mergeList_
145 [
146 blockCell[cellPointi]
147 + blockOffsets_[blocki]
148 ];
149 }
150
151 // Construct collapsed cell and add to list
152 cells_[celli].reset(hex, cellPoints, true);
153 ++celli;
154 }
155 }
156}
157
158
159Foam::faceList Foam::blockMesh::createPatchFaces
160(
161 const polyPatch& patchTopologyFaces
162) const
163{
164 const blockList& blocks = *this;
165
166 labelList blockLabels = patchTopologyFaces.polyPatch::faceCells();
167
168 label nFaces = 0;
169
170 forAll(patchTopologyFaces, patchTopologyFaceLabel)
171 {
172 const label blocki = blockLabels[patchTopologyFaceLabel];
173
174 faceList blockFaces = blocks[blocki].blockShape().faces();
175
176 forAll(blockFaces, blockFaceLabel)
177 {
178 if
179 (
180 blockFaces[blockFaceLabel]
181 == patchTopologyFaces[patchTopologyFaceLabel]
182 )
183 {
184 nFaces +=
185 blocks[blocki].boundaryPatches()[blockFaceLabel].size();
186 }
187 }
188 }
189
190
191 faceList patchFaces(nFaces);
192 face quad(4);
193 label faceLabel = 0;
194
195 forAll(patchTopologyFaces, patchTopologyFaceLabel)
196 {
197 const label blocki = blockLabels[patchTopologyFaceLabel];
198
199 faceList blockFaces = blocks[blocki].blockShape().faces();
200
201 forAll(blockFaces, blockFaceLabel)
202 {
203 if
204 (
205 blockFaces[blockFaceLabel]
206 == patchTopologyFaces[patchTopologyFaceLabel]
207 )
208 {
209 const List<FixedList<label, 4>>& blockPatchFaces =
210 blocks[blocki].boundaryPatches()[blockFaceLabel];
211
212 forAll(blockPatchFaces, blockFaceLabel)
213 {
214 // Lookup the face points
215 // and collapse duplicate point labels
216
217 quad[0] =
218 mergeList_
219 [
220 blockPatchFaces[blockFaceLabel][0]
221 + blockOffsets_[blocki]
222 ];
223
224 label nUnique = 1;
225
226 for
227 (
228 label facePointLabel = 1;
229 facePointLabel < 4;
230 facePointLabel++
231 )
232 {
233 quad[nUnique] =
234 mergeList_
235 [
236 blockPatchFaces[blockFaceLabel][facePointLabel]
237 + blockOffsets_[blocki]
238 ];
239
240 if (quad[nUnique] != quad[nUnique-1])
241 {
242 nUnique++;
243 }
244 }
245
246 if (quad[nUnique-1] == quad[0])
247 {
248 nUnique--;
249 }
250
251 if (nUnique == 4)
252 {
253 patchFaces[faceLabel++] = quad;
254 }
255 else if (nUnique == 3)
256 {
257 patchFaces[faceLabel++] = face
258 (
259 labelList::subList(quad, 3)
260 );
261 }
262 // else the face has collapsed to an edge or point
263 }
264 }
265 }
266 }
267
268 patchFaces.resize(faceLabel);
269
270 return patchFaces;
271}
272
273
274void Foam::blockMesh::createPatches() const
275{
276 const polyPatchList& topoPatches = topology().boundaryMesh();
277
278 if (verbose_)
279 {
280 Info<< "Creating patches" << endl;
281 }
282
283 patches_.resize(topoPatches.size());
284
285 forAll(topoPatches, patchi)
286 {
287 patches_[patchi] = createPatchFaces(topoPatches[patchi]);
288 }
289}
290
291
292// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
293
295{
296 if (!topologyPtr_)
297 {
299 << "topology not allocated"
300 << exit(FatalError);
301 }
302
303 return *topologyPtr_;
304}
305
306
308Foam::blockMesh::topology(bool applyTransform) const
309{
310 const polyMesh& origTopo = topology();
311
312 if (applyTransform && hasPointTransforms())
313 {
314 // Copy mesh components
315
317 newIO.registerObject(false);
318
319 pointField newPoints(origTopo.points());
320 inplacePointTransforms(newPoints);
321
322 auto ttopoMesh = refPtr<polyMesh>::New
323 (
324 newIO,
325 std::move(newPoints),
326 faceList(origTopo.faces()),
327 labelList(origTopo.faceOwner()),
328 labelList(origTopo.faceNeighbour())
329 );
330 auto& topoMesh = ttopoMesh.ref();
331
332 // Patches
333 const polyBoundaryMesh& pbmOld = origTopo.boundaryMesh();
334 const polyBoundaryMesh& pbmNew = topoMesh.boundaryMesh();
335
336 polyPatchList newPatches(pbmOld.size());
337
338 forAll(pbmOld, patchi)
339 {
340 newPatches.set(patchi, pbmOld[patchi].clone(pbmNew));
341 }
342
343 topoMesh.addPatches(newPatches);
344
345 return ttopoMesh;
346 }
347 else
348 {
349 return origTopo;
350 }
351}
352
353
356{
357 const blockMesh& blkMesh = *this;
358
359 if (verbose_)
360 {
361 Info<< nl << "Creating polyMesh from blockMesh" << endl;
362 }
363
365 (
366 io,
367 pointField(blkMesh.points()), // Use a copy of the block points
368 blkMesh.cells(),
369 blkMesh.patches(),
370 blkMesh.patchNames(),
371 blkMesh.patchDicts(),
372 "defaultFaces", // Default patch name
373 emptyPolyPatch::typeName // Default patch type
374 );
375
376
377 // Set any cellZones
378 const label nZones = blkMesh.numZonedBlocks();
379
380 if (nZones)
381 {
382 polyMesh& pmesh = *meshPtr;
383
384 if (verbose_)
385 {
386 Info<< "Adding cell zones" << endl;
387 }
388
389 // Map from zoneName to cellZone index
390 HashTable<label> zoneMap(2*nZones);
391
392 // Cells per zone
394
395 // Running cell counter
396 label celli = 0;
397
398 // Largest zone so far
399 label freeZonei = 0;
400
401 for (const block& b : blkMesh)
402 {
403 const word& zoneName = b.zoneName();
404 const label nCellsInBlock = b.cells().size();
405
406 if (zoneName.size())
407 {
408 const auto iter = zoneMap.cfind(zoneName);
409
410 label zonei = freeZonei;
411
412 if (iter.found())
413 {
414 zonei = *iter;
415 }
416 else
417 {
418 zoneMap.insert(zoneName, zonei);
419 ++freeZonei;
420
421 if (verbose_)
422 {
423 Info<< " " << zonei << '\t' << zoneName << endl;
424 }
425 }
426
427
428 // Fill with cell ids
429
430 zoneCells[zonei].reserve
431 (
432 zoneCells[zonei].size() + nCellsInBlock
433 );
434
435 const label endOfFill = celli + nCellsInBlock;
436
437 for (; celli < endOfFill; ++celli)
438 {
439 zoneCells[zonei].append(celli);
440 }
441 }
442 else
443 {
444 celli += nCellsInBlock;
445 }
446 }
447
448 List<cellZone*> cz(zoneMap.size());
449 forAllConstIters(zoneMap, iter)
450 {
451 const word& zoneName = iter.key();
452 const label zonei = iter.val();
453
454 cz[zonei] = new cellZone
455 (
456 zoneName,
457 zoneCells[zonei].shrink(),
458 zonei,
459 pmesh.cellZones()
460 );
461 }
462
463 pmesh.pointZones().clear();
464 pmesh.faceZones().clear();
465 pmesh.cellZones().clear();
467 }
468
469
470 // Merge patch pairs, cyclic must be done elsewhere
471 // - requires libdynamicMesh
472
473 return meshPtr;
474}
475
476
477// ************************************************************************* //
const polyMesh & topoMesh
Definition: blockMeshOBJ.H:29
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition: HashTableI.H:141
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
bool registerObject() const noexcept
Should object created with this IOobject be registered?
Definition: IOobjectI.H:107
SubList< label > subList
Declare type of subList.
Definition: List.H:111
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
void clear()
Clear the zones.
Definition: ZoneMesh.C:730
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A multi-block mesh generator.
Definition: blockMesh.H:168
const faceListList & patches() const
Return the patch face lists.
Definition: blockMesh.C:388
wordList patchNames() const
Return the patch names.
Definition: blockMesh.C:399
label numZonedBlocks() const
Number of blocks with specified zones.
Definition: blockMesh.C:417
bool inplacePointTransforms(pointField &pts) const
Apply coordinate transforms and scaling.
Definition: blockMesh.C:441
const cellShapeList & cells() const
Return cell shapes list.
Definition: blockMesh.C:377
const pointField & points() const
Definition: blockMesh.C:366
PtrList< dictionary > patchDicts() const
Patch information from the topology mesh.
Definition: blockMesh.C:349
const polyMesh & topology() const
PtrList< block > blockList
The list of blocks is stored as a PtrList.
Definition: blockMesh.H:174
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:61
A subset of mesh cells.
Definition: cellZone.H:65
reference ref() const
A reference to the entry (Error if not found)
Definition: dictionary.H:225
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:498
void addPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:975
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:504
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1127
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:1013
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:492
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
A class for managing references or pointers (no reference counting)
Definition: refPtr.H:58
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition: word.H:68
const cellModel & hex
Foam::autoPtr< Foam::dynamicFvMesh > meshPtr
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
const cellShapeList & cells
label nZones
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition: polyPatch.H:63
List< label > labelList
A List of labels.
Definition: List.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
PtrList< block > blockList
A PtrList of blocks.
Definition: blockList.H:47
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
SubList< label > labelSubList
A SubList of labels.
Definition: SubList.H:59
error FatalError
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
static const char *const typeName
The type name used in ensight case files.