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 -------------------------------------------------------------------------------
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  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 
121 void 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 
159 Foam::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 
274 void 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 
308 Foam::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  PtrList<polyPatch> 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
393  List<DynamicList<label>> zoneCells(nZones);
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();
466  pmesh.addZones(List<pointZone*>(), List<faceZone*>(), cz);
467  }
468 
469 
470  // Merge patch pairs, cyclic must be done elsewhere
471  // - requires libdynamicMesh
472 
473  return meshPtr;
474 }
475 
476 
477 // ************************************************************************* //
nZones
label nZones
Definition: interpolatedFaces.H:24
Foam::refPtr::New
static refPtr< T > New(Args &&... args)
Construct refPtr of T with forwarding arguments.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::polyMesh::addPatches
void addPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:961
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:67
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::ZoneMesh::clear
void clear()
Clear the zones.
Definition: ZoneMesh.C:724
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::blockMesh::patchDicts
PtrList< dictionary > patchDicts() const
Patch information from the topology mesh.
Definition: blockMesh.C:369
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
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:355
Foam::blockMesh::patches
const faceListList & patches() const
Return the patch face lists.
Definition: blockMesh.C:408
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::IOobject::registerObject
bool registerObject() const noexcept
Should object created with this IOobject be registered?
Definition: IOobjectI.H:107
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::cellModel::HEX
hex
Definition: cellModel.H:81
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::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
blockMesh.H
Foam::blockMesh::topology
const polyMesh & topology() const
Definition: blockMeshCreate.C:294
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
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
topoMesh
const polyMesh & topoMesh
Definition: blockMeshOBJ.H:29
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
cellModel.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::labelSubList
SubList< label > labelSubList
A SubList of labels.
Definition: SubList.H:59
Foam::polyPatchList
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:47
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:492
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:486
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
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::blockMesh::inplacePointTransforms
bool inplacePointTransforms(pointField &pts) const
Apply coordinate transforms and scaling.
Definition: blockMesh.C:461
Foam::List::subList
SubList< T > subList
Declare type of subList.
Definition: List.H:112
Foam::FatalError
error FatalError
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
emptyPolyPatch.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::hex
IOstream & hex(IOstream &io)
Definition: IOstream.H:446
Foam::HashTable< label >
Foam::blockMesh::blockList
PtrList< block > blockList
The list of blocks is stored as a PtrList.
Definition: blockMesh.H:174
Foam::autoPtr< Foam::polyMesh >
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::blockMesh::numZonedBlocks
label numZonedBlocks() const
Number of blocks with specified zones.
Definition: blockMesh.C:437
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:404
Foam::blockMesh::points
const pointField & points() const
Definition: blockMesh.C:386
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 the patch names.
Definition: blockMesh.C:419
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::polyMesh::pointZones
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:480
Foam::blockList
PtrList< block > blockList
A PtrList of blocks.
Definition: blockList.H:47
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::blockMesh::cells
const cellShapeList & cells() const
Return cell shapes list.
Definition: blockMesh.C:397
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::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::refPtr
A class for managing references or pointers (no reference counting)
Definition: PtrList.H:60
Foam::blockMesh
A multi-block mesh generator.
Definition: blockMesh.H:165