blockMeshMergeTopological.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) 2015-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 
31 // * * * * * * * * * * * * * * * Private Functions * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // faces
37 // 6
38 // (
39 // 4(0 4 7 3) // x-min
40 // 4(1 2 6 5) // x-max
41 // 4(0 1 5 4) // y-min
42 // 4(3 7 6 2) // y-max
43 // 4(0 3 2 1) // z-min
44 // 4(4 5 6 7) // z-max
45 // );
46 
47 // Face-edge directions
48 static const int faceEdgeDirs[6][4] =
49 {
50  {2, 1, -2, -1},
51  {1, 2, -1, -2},
52  {1, 2, -1, -2},
53  {2, 1, -2, -1},
54  {2, 1, -2, -1},
55  {1, 2, -1, -2}
56 };
57 
58 // The face-face-rotation direction correspondence map
59 static Pair<int> faceFaceRotMap[6][6][4];
60 
61 // Generate the face-face-rotation direction correspondence map
63 {
64  for(int facePi=0; facePi<6; facePi++)
65  {
66  for(int faceNi=0; faceNi<6; faceNi++)
67  {
68  for(int rot=0; rot<4; rot++)
69  {
70  Pair<int>& map = faceFaceRotMap[facePi][faceNi][rot];
71 
72  for(int Pp=0; Pp<2; Pp++)
73  {
74  int Pdir = faceEdgeDirs[facePi][Pp];
75  int Np = (3 - Pp + rot)%4;
76  int Ndir = faceEdgeDirs[faceNi][Np];
77  map[Pdir-1] = -Ndir;
78  }
79 
80  // Handle sign change due to match-face transpose
81  if (mag(map[0]) == 2 && map[0]*map[1] < 0)
82  {
83  map[0] = -map[0];
84  map[1] = -map[1];
85  }
86  }
87  }
88  }
89 }
90 
91 
92 // Return the direction map for the merge-faces
93 Pair<int> faceMap
94 (
95  const label facePi,
96  const face& faceP,
97  const label faceNi,
98  const face& faceN
99 )
100 {
101  // Search for the point on faceN corresponding to the 0-point on faceP
102  for(int rot=0; rot<4; rot++)
103  {
104  if (faceN[rot] == faceP[0])
105  {
106  return faceFaceRotMap[facePi][faceNi][rot];
107  }
108  }
109 
111  << "Cannot find point correspondence for faces "
112  << faceP << " and " << faceN
113  << exit(FatalError);
114 
115  return Pair<int>(0, 0);
116 }
117 
118 
119 // Set the block and face indices for all the merge faces
121 (
122  const cellList& topoCells,
123  const faceList::subList& topoInternalFaces,
124  const labelList& topoFaceCell,
125  List<Pair<label>>& mergeBlock
126 )
127 {
128  forAll(topoInternalFaces, topoFacei)
129  {
130  label topoPi = topoFaceCell[topoFacei];
131  const labelList& topoPfaces = topoCells[topoPi];
132 
133  bool foundFace = false;
134  label topoPfacei;
135  for
136  (
137  topoPfacei = 0;
138  topoPfacei < topoPfaces.size();
139  topoPfacei++
140  )
141  {
142  if (topoPfaces[topoPfacei] == topoFacei)
143  {
144  foundFace = true;
145  break;
146  }
147  }
148 
149  if (!foundFace)
150  {
152  << "Cannot find merge face for block " << topoPi
153  << exit(FatalError);
154  }
155 
156  mergeBlock[topoFacei].first() = topoPi;
157  mergeBlock[topoFacei].second() = topoPfacei;
158  }
159 }
160 
161 
162 // Return the number of divisions in each direction for the face
163 Pair<label> faceNij(const label facei, const block& block)
164 {
165  Pair<label> fnij;
166 
167  int i = facei/2;
168 
169  if (i == 0)
170  {
171  fnij.first() = block.density().y() + 1;
172  fnij.second() = block.density().z() + 1;
173  }
174  else if (i == 1)
175  {
176  fnij.first() = block.density().x() + 1;
177  fnij.second() = block.density().z() + 1;
178  }
179  else if (i == 2)
180  {
181  fnij.first() = block.density().x() + 1;
182  fnij.second() = block.density().y() + 1;
183  }
184 
185  return fnij;
186 }
187 
188 
189 // Sign the index corresponding to the map
190 inline label signIndex(const int map, const label i)
191 {
192  return map < 0 ? -i-1 : i;
193 }
194 
195 
196 // Reverse a signed index with the number of divisions
197 inline label unsignIndex(const label i, const label ni)
198 {
199  return i >= 0 ? i : ni + i + 1;
200 }
201 
202 
203 // Return the mapped index
204 inline label mapij(const int map, const label i, const label j)
205 {
206  return signIndex(map, mag(map) == 1 ? i : j);
207 }
208 
209 
210 // Return the face point index
211 inline label facePoint
212 (
213  const int facei,
214  const block& block,
215  const label i,
216  const label j
217 )
218 {
219  switch (facei)
220  {
221  case 0:
222  return block.pointLabel(0, i, j);
223  case 1:
224  return block.pointLabel(block.density().x(), i, j);
225  case 2:
226  return block.pointLabel(i, 0, j);
227  case 3:
228  return block.pointLabel(i, block.density().y(), j);
229  case 4:
230  return block.pointLabel(i, j, 0);
231  case 5:
232  return block.pointLabel(i, j, block.density().z());
233  default:
234  return -1;
235  }
236 }
237 
238 
239 // Return the neighbour face point from the signed indices
240 inline label facePointN
241 (
242  const block& block,
243  const label i,
244  const label j,
245  const label k
246 )
247 {
248  return block.pointLabel
249  (
250  unsignIndex(i, block.density().x()),
251  unsignIndex(j, block.density().y()),
252  unsignIndex(k, block.density().z())
253  );
254 }
255 
256 
257 // Return the neighbour face point from the mapped indices
258 inline label facePointN
259 (
260  const int facei,
261  const block& block,
262  const label i,
263  const label j
264 )
265 {
266  switch (facei)
267  {
268  case 0:
269  return facePointN(block, 0, i, j);
270  case 1:
271  return facePointN(block, block.density().x(), i, j);
272  case 2:
273  return facePointN(block, i, 0, j);
274  case 3:
275  return facePointN(block, i, block.density().y(), j);
276  case 4:
277  return facePointN(block, i, j, 0);
278  case 5:
279  return facePointN(block, i, j, block.density().z());
280  default:
281  return -1;
282  }
283 }
284 
285 
286 // Return the neighbour face point using the map
287 inline label facePointN
288 (
289  const int facei,
290  const Pair<int>& fmap,
291  const block& block,
292  const label i,
293  const label j
294 )
295 {
296  return facePointN(facei, block, mapij(fmap[0], i, j), mapij(fmap[1], i, j));
297 }
298 
299 } // End namespace Foam
300 
301 
302 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
303 
304 void Foam::blockMesh::calcTopologicalMerge()
305 {
306  // Generate the static face-face map
308 
309  const blockList& blocks = *this;
310 
311  if (verbose_)
312  {
313  Info<< "Creating block offsets" << endl;
314  }
315 
316  blockOffsets_.setSize(blocks.size());
317 
318  nPoints_ = 0;
319  nCells_ = 0;
320 
321  forAll(blocks, blocki)
322  {
323  blockOffsets_[blocki] = nPoints_;
324 
325  nPoints_ += blocks[blocki].nPoints();
326  nCells_ += blocks[blocki].nCells();
327  }
328 
329  if (verbose_)
330  {
331  Info<< "Creating merge list (topological search).." << flush;
332  }
333 
334  // Size merge list and initialize to -1
335  mergeList_.setSize(nPoints_, -1);
336 
337  // Block mesh topology
338  const pointField& topoPoints = topology().points();
339  const cellList& topoCells = topology().cells();
340  const faceList& topoFaces = topology().faces();
341  const labelList& topoFaceOwn = topology().faceOwner();
342  const labelList& topoFaceNei = topology().faceNeighbour();
343 
344  // Topological merging is only necessary for internal block faces
345  // Note edge and face collapse may apply to boundary faces
346  // but is not yet supported in the "fast" algorithm
347  const faceList::subList topoInternalFaces
348  (
349  topoFaces,
350  topology().nInternalFaces()
351  );
352 
353  List<Pair<label>> mergeBlockP(topoInternalFaces.size());
355  (
356  topoCells,
357  topoInternalFaces,
358  topoFaceOwn,
359  mergeBlockP
360  );
361 
362  List<Pair<label>> mergeBlockN(topoInternalFaces.size());
364  (
365  topoCells,
366  topoInternalFaces,
367  topoFaceNei,
368  mergeBlockN
369  );
370 
371  DebugInfo << endl;
372 
373  forAll(topoInternalFaces, topoFacei)
374  {
375  DebugInfo << "Processing face " << topoFacei << endl;
376 
377  label blockPi = mergeBlockP[topoFacei].first();
378  label blockPfacei = mergeBlockP[topoFacei].second();
379 
380  label blockNi = mergeBlockN[topoFacei].first();
381  label blockNfacei = mergeBlockN[topoFacei].second();
382 
383  Pair<int> fmap
384  (
385  faceMap
386  (
387  blockPfacei,
388  blocks[blockPi].blockShape().faces()[blockPfacei],
389  blockNfacei,
390  blocks[blockNi].blockShape().faces()[blockNfacei]
391  )
392  );
393 
394  DebugInfo
395  << " Face map for faces "
396  << blocks[blockPi].blockShape().faces()[blockPfacei] << " "
397  << blocks[blockNi].blockShape().faces()[blockNfacei] << ": "
398  << fmap << endl;
399 
400  const pointField& blockPpoints = blocks[blockPi].points();
401  const pointField& blockNpoints = blocks[blockNi].points();
402 
403  Pair<label> Pnij(faceNij(blockPfacei, blocks[blockPi]));
404 
405  // Check block subdivision correspondence
406  {
407  Pair<label> Nnij(faceNij(blockNfacei, blocks[blockNi]));
408  Pair<label> NPnij;
409  NPnij[0] = Nnij[mag(fmap[0]) - 1];
410  NPnij[1] = Nnij[mag(fmap[1]) - 1];
411 
412  if (Pnij != NPnij)
413  {
415  << "Sub-division mismatch between face "
416  << blockPfacei << " of block " << blockPi << Pnij
417  << " and face "
418  << blockNfacei << " of block " << blockNi << Nnij
419  << exit(FatalError);
420  }
421  }
422 
423  // Calculate a suitable test distance from the bounding box of the face.
424  // Note this is used only as a sanity check and for diagnostics in
425  // case there is a grading inconsistency.
426  const boundBox bb(topoCells[blockPi].points(topoFaces, topoPoints));
427  const scalar testSqrDist = magSqr(1e-6*bb.span());
428 
429  // Accumulate the maximum merge distance for diagnostics
430  scalar maxSqrDist = 0;
431 
432  for (label j=0; j<Pnij.second(); j++)
433  {
434  for (label i=0; i<Pnij.first(); i++)
435  {
436  label blockPpointi =
437  facePoint(blockPfacei, blocks[blockPi], i, j);
438 
439  label blockNpointi =
440  facePointN(blockNfacei, fmap, blocks[blockNi], i, j);
441 
442  scalar sqrDist
443  (
444  magSqr
445  (
446  blockPpoints[blockPpointi]
447  - blockNpoints[blockNpointi]
448  )
449  );
450 
451  if (sqrDist > testSqrDist)
452  {
454  << "Point merge failure between face "
455  << blockPfacei << " of block " << blockPi
456  << " and face "
457  << blockNfacei << " of block " << blockNi
458  << endl
459  << " Points: " << blockPpoints[blockPpointi]
460  << " " << blockNpoints[blockNpointi]
461  << endl
462  << " This may be due to inconsistent grading."
463  << exit(FatalError);
464  }
465 
466  maxSqrDist = max(maxSqrDist, sqrDist);
467 
468  label Ppointi = blockPpointi + blockOffsets_[blockPi];
469  label Npointi = blockNpointi + blockOffsets_[blockNi];
470 
471  label minPNi = min(Ppointi, Npointi);
472 
473  if (mergeList_[Ppointi] != -1)
474  {
475  minPNi = min(minPNi, mergeList_[Ppointi]);
476  }
477 
478  if (mergeList_[Npointi] != -1)
479  {
480  minPNi = min(minPNi, mergeList_[Npointi]);
481  }
482 
483  mergeList_[Ppointi] = mergeList_[Npointi] = minPNi;
484  }
485  }
486 
487  DebugInfo
488  << " Max distance between merge points: "
489  << sqrt(maxSqrDist) << endl;
490  }
491 
492 
493  bool changedPointMerge = false;
494  label nPasses = 0;
495 
496  do
497  {
498  changedPointMerge = false;
499  nPasses++;
500 
501  forAll(topoInternalFaces, topoFacei)
502  {
503  label blockPi = mergeBlockP[topoFacei].first();
504  label blockPfacei = mergeBlockP[topoFacei].second();
505 
506  label blockNi = mergeBlockN[topoFacei].first();
507  label blockNfacei = mergeBlockN[topoFacei].second();
508 
509  Pair<int> fmap
510  (
511  faceMap
512  (
513  blockPfacei,
514  blocks[blockPi].blockShape().faces()[blockPfacei],
515  blockNfacei,
516  blocks[blockNi].blockShape().faces()[blockNfacei]
517  )
518  );
519 
520  Pair<label> Pnij(faceNij(blockPfacei, blocks[blockPi]));
521 
522  for (label j=0; j<Pnij.second(); j++)
523  {
524  for (label i=0; i<Pnij.first(); i++)
525  {
526  label blockPpointi =
527  facePoint(blockPfacei, blocks[blockPi], i, j);
528 
529  label blockNpointi =
530  facePointN(blockNfacei, fmap, blocks[blockNi], i, j);
531 
532  label Ppointi =
533  blockPpointi + blockOffsets_[blockPi];
534 
535  label Npointi =
536  blockNpointi + blockOffsets_[blockNi];
537 
538  if (mergeList_[Ppointi] != mergeList_[Npointi])
539  {
540  changedPointMerge = true;
541 
542  mergeList_[Ppointi]
543  = mergeList_[Npointi]
544  = min(mergeList_[Ppointi], mergeList_[Npointi]);
545  }
546  }
547  }
548  }
549 
550  if (verbose_)
551  {
552  Info<< '.' << flush;
553  }
554 
555  if (nPasses > 100)
556  {
558  << "Point merging failed after 100 passes."
559  << exit(FatalError);
560  }
561 
562  } while (changedPointMerge);
563 
564  if (verbose_)
565  {
566  Info<< endl;
567  }
568 
569  // Sort merge list and count number of unique points
570  label nUniqPoints = 0;
571 
572  forAll(mergeList_, pointi)
573  {
574  if (mergeList_[pointi] > pointi)
575  {
577  << "Merge list contains point index out of range"
578  << exit(FatalError);
579  }
580 
581  if (mergeList_[pointi] == -1 || mergeList_[pointi] == pointi)
582  {
583  mergeList_[pointi] = nUniqPoints++;
584  }
585  else
586  {
587  mergeList_[pointi] = mergeList_[mergeList_[pointi]];
588  }
589  }
590 
591  // Correct number of points in mesh
592  nPoints_ = nUniqPoints;
593 }
594 
595 
596 // ************************************************************************* //
Foam::Pair::second
const T & second() const noexcept
Return second element, which is also the last element.
Definition: PairI.H:122
Foam::unsignIndex
label unsignIndex(const label i, const label ni)
Definition: blockMeshMergeTopological.C:197
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::faceNij
Pair< label > faceNij(const label facei, const block &block)
Definition: blockMeshMergeTopological.C:163
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
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
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::setBlockFaceCorrespondence
void setBlockFaceCorrespondence(const cellList &topoCells, const faceList::subList &topoInternalFaces, const labelList &topoFaceCell, List< Pair< label >> &mergeBlock)
Definition: blockMeshMergeTopological.C:121
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
blockMesh.H
Foam::mapij
label mapij(const int map, const label i, const label j)
Definition: blockMeshMergeTopological.C:204
Foam::blockMesh::faces
const blockFaceList & faces() const
Return the curved faces.
Definition: blockMesh.H:351
Foam::blockMesh::topology
const polyMesh & topology() const
Return the blockMesh topology as a polyMesh.
Definition: blockMesh.C:140
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::facePoint
label facePoint(const int facei, const block &block, const label i, const label j)
Definition: blockMeshMergeTopological.C:212
Foam::flush
Ostream & flush(Ostream &os)
Flush stream.
Definition: Ostream.H:342
Foam::ijkMesh::pointLabel
label pointLabel(const label i, const label j, const label k) const
The linear point index for an i-j-k position.
Definition: ijkMeshI.H:183
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::cellList
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::List< face >::subList
SubList< face > subList
Declare type of subList.
Definition: List.H:114
Foam::FatalError
error FatalError
topoCells
const vtk::vtuCells topoCells(topoMesh, writeOpts)
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::faceFaceRotMap
static Pair< int > faceFaceRotMap[6][6][4]
Definition: blockMeshMergeTopological.C:59
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::blockMesh::blockList
PtrList< block > blockList
The list of blocks is stored as a PtrList.
Definition: blockMesh.H:157
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:381
Foam::facePointN
label facePointN(const block &block, const label i, const label j, const label k)
Definition: blockMeshMergeTopological.C:241
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:359
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::blockMesh::points
const pointField & points() const
The points for the entire mesh.
Definition: blockMesh.C:176
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::signIndex
label signIndex(const int map, const label i)
Definition: blockMeshMergeTopological.C:190
Foam::faceEdgeDirs
static const int faceEdgeDirs[6][4]
Definition: blockMeshMergeTopological.C:48
Foam::List< cell >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::blockDescriptor::density
const labelVector & density() const
Return the mesh density (number of cells) in the i,j,k directions.
Definition: blockDescriptorI.H:49
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::genFaceFaceRotMap
void genFaceFaceRotMap()
Definition: blockMeshMergeTopological.C:62