meshDualiser.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 "meshDualiser.H"
27 #include "meshTools.H"
28 #include "polyMesh.H"
29 #include "polyTopoChange.H"
30 #include "polyTopoChangeMap.H"
31 #include "edgeFaceCirculator.H"
32 #include "OFstream.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(meshDualiser, 0);
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 // Dump state so far.
45 void Foam::meshDualiser::dumpPolyTopoChange
46 (
47  const polyTopoChange& meshMod,
48  const fileName& prefix
49 )
50 {
51  OFstream str1(prefix + "Faces.obj");
52  OFstream str2(prefix + "Edges.obj");
53 
54  Info<< "Dumping current polyTopoChange. Faces to " << str1.name()
55  << " , points and edges to " << str2.name() << endl;
56 
57  const DynamicList<point>& points = meshMod.points();
58 
59  forAll(points, pointi)
60  {
61  meshTools::writeOBJ(str1, points[pointi]);
62  meshTools::writeOBJ(str2, points[pointi]);
63  }
64 
65  const DynamicList<face>& faces = meshMod.faces();
66 
67  forAll(faces, facei)
68  {
69  const face& f = faces[facei];
70 
71  str1<< 'f';
72  forAll(f, fp)
73  {
74  str1<< ' ' << f[fp]+1;
75  }
76  str1<< nl;
77 
78  str2<< 'l';
79  forAll(f, fp)
80  {
81  str2<< ' ' << f[fp]+1;
82  }
83  str2<< ' ' << f[0]+1 << nl;
84  }
85 }
86 
87 
88 Foam::label Foam::meshDualiser::findDualCell
89 (
90  const label celli,
91  const label pointi
92 ) const
93 {
94  const labelList& dualCells = pointToDualCells_[pointi];
95 
96  if (dualCells.size() == 1)
97  {
98  return dualCells[0];
99  }
100  else
101  {
102  label index = findIndex(mesh_.pointCells()[pointi], celli);
103 
104  return dualCells[index];
105  }
106 }
107 
108 
109 void Foam::meshDualiser::generateDualBoundaryEdges
110 (
111  const PackedBoolList& isBoundaryEdge,
112  const label pointi,
113  polyTopoChange& meshMod
114 )
115 {
116  const labelList& pEdges = mesh_.pointEdges()[pointi];
117 
118  forAll(pEdges, pEdgeI)
119  {
120  label edgeI = pEdges[pEdgeI];
121 
122  if (edgeToDualPoint_[edgeI] == -1 && isBoundaryEdge.get(edgeI) == 1)
123  {
124  const edge& e = mesh_.edges()[edgeI];
125 
126  edgeToDualPoint_[edgeI] = meshMod.addPoint
127  (
128  e.centre(mesh_.points()),
129  pointi, // masterPoint
130  -1, // zoneID
131  true // inCell
132  );
133  }
134  }
135 }
136 
137 
138 // Return true if point on face has same dual cells on both owner and neighbour
139 // sides.
140 bool Foam::meshDualiser::sameDualCell
141 (
142  const label facei,
143  const label pointi
144 ) const
145 {
146  if (!mesh_.isInternalFace(facei))
147  {
149  << "face:" << facei << " is not internal face."
150  << abort(FatalError);
151  }
152 
153  label own = mesh_.faceOwner()[facei];
154  label nei = mesh_.faceNeighbour()[facei];
155 
156  return findDualCell(own, pointi) == findDualCell(nei, pointi);
157 }
158 
159 
160 Foam::label Foam::meshDualiser::addInternalFace
161 (
162  const label masterPointi,
163  const label masterEdgeI,
164  const label masterFacei,
165 
166  const bool edgeOrder,
167  const label dualCell0,
168  const label dualCell1,
169  const DynamicList<label>& verts,
170  polyTopoChange& meshMod
171 ) const
172 {
173  face newFace(verts);
174 
175  if (edgeOrder != (dualCell0 < dualCell1))
176  {
177  reverse(newFace);
178  }
179 
180  label zoneID = -1;
181  bool zoneFlip = false;
182  if (masterFacei != -1)
183  {
184  zoneID = mesh_.faceZones().whichZone(masterFacei);
185 
186  if (zoneID != -1)
187  {
188  const faceZone& fZone = mesh_.faceZones()[zoneID];
189 
190  zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
191  }
192  }
193 
194  label dualFacei;
195 
196  if (dualCell0 < dualCell1)
197  {
198  dualFacei = meshMod.addFace
199  (
200  newFace,
201  dualCell0, // own
202  dualCell1, // nei
203  masterPointi, // masterPointID
204  masterEdgeI, // masterEdgeID
205  masterFacei, // masterFaceID
206  false, // flipFaceFlux
207  -1, // patchID
208  zoneID, // zoneID
209  zoneFlip // zoneFlip
210  );
211 
212  // pointField dualPoints(meshMod.points());
213  // vector n(newFace.normal(dualPoints));
214  // n /= mag(n);
215  // Pout<< "Generated internal dualFace:" << dualFacei
216  // << " verts:" << newFace
217  // << " points:" << UIndirectList<point>(meshMod.points(), newFace)()
218  // << " n:" << n
219  // << " between dualowner:" << dualCell0
220  // << " dualneighbour:" << dualCell1
221  // << endl;
222  }
223  else
224  {
225  dualFacei = meshMod.addFace
226  (
227  newFace,
228  dualCell1, // own
229  dualCell0, // nei
230  masterPointi, // masterPointID
231  masterEdgeI, // masterEdgeID
232  masterFacei, // masterFaceID
233  false, // flipFaceFlux
234  -1, // patchID
235  zoneID, // zoneID
236  zoneFlip // zoneFlip
237  );
238 
239  // pointField dualPoints(meshMod.points());
240  // vector n(newFace.normal(dualPoints));
241  // n /= mag(n);
242  // Pout<< "Generated internal dualFace:" << dualFacei
243  // << " verts:" << newFace
244  // << " points:" << UIndirectList<point>(meshMod.points(), newFace)()
245  // << " n:" << n
246  // << " between dualowner:" << dualCell1
247  // << " dualneighbour:" << dualCell0
248  // << endl;
249  }
250  return dualFacei;
251 }
252 
253 
254 Foam::label Foam::meshDualiser::addBoundaryFace
255 (
256  const label masterPointi,
257  const label masterEdgeI,
258  const label masterFacei,
259 
260  const label dualCelli,
261  const label patchi,
262  const DynamicList<label>& verts,
263  polyTopoChange& meshMod
264 ) const
265 {
266  face newFace(verts);
267 
268  label zoneID = -1;
269  bool zoneFlip = false;
270  if (masterFacei != -1)
271  {
272  zoneID = mesh_.faceZones().whichZone(masterFacei);
273 
274  if (zoneID != -1)
275  {
276  const faceZone& fZone = mesh_.faceZones()[zoneID];
277 
278  zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
279  }
280  }
281 
282  label dualFacei = meshMod.addFace
283  (
284  newFace,
285  dualCelli, // own
286  -1, // nei
287  masterPointi, // masterPointID
288  masterEdgeI, // masterEdgeID
289  masterFacei, // masterFaceID
290  false, // flipFaceFlux
291  patchi, // patchID
292  zoneID, // zoneID
293  zoneFlip // zoneFlip
294  );
295 
296  // pointField dualPoints(meshMod.points());
297  // vector n(newFace.normal(dualPoints));
298  // n /= mag(n);
299  // Pout<< "Generated boundary dualFace:" << dualFacei
300  // << " verts:" << newFace
301  // << " points:" << UIndirectList<point>(meshMod.points(), newFace)()
302  // << " n:" << n
303  // << " on dualowner:" << dualCelli
304  // << endl;
305  return dualFacei;
306 }
307 
308 
309 // Walks around edgeI.
310 // splitFace=true : creates multiple faces
311 // splitFace=false: creates single face if same dual cells on both sides,
312 // multiple faces otherwise.
313 void Foam::meshDualiser::createFacesAroundEdge
314 (
315  const bool splitFace,
316  const PackedBoolList& isBoundaryEdge,
317  const label edgeI,
318  const label startFacei,
319  polyTopoChange& meshMod,
320  boolList& doneEFaces
321 ) const
322 {
323  const edge& e = mesh_.edges()[edgeI];
324  const labelList& eFaces = mesh_.edgeFaces()[edgeI];
325 
327  (
328  mesh_.faces()[startFacei],
329  e[0],
330  e[1]
331  );
332 
333  edgeFaceCirculator ie
334  (
335  mesh_,
336  startFacei, // face
337  true, // ownerSide
338  fp, // fp
339  isBoundaryEdge.get(edgeI) == 1 // isBoundaryEdge
340  );
341  ie.setCanonical();
342 
343  bool edgeOrder = ie.sameOrder(e[0], e[1]);
344  label startFaceLabel = ie.faceLabel();
345 
346  // Pout<< "At edge:" << edgeI << " verts:" << e
347  // << " points:" << mesh_.points()[e[0]] << mesh_.points()[e[1]]
348  // << " started walking at face:" << ie.faceLabel()
349  // << " verts:" << mesh_.faces()[ie.faceLabel()]
350  // << " edgeOrder:" << edgeOrder
351  // << " in direction of cell:" << ie.cellLabel()
352  // << endl;
353 
354  // Walk and collect face.
355  DynamicList<label> verts(100);
356 
357  if (edgeToDualPoint_[edgeI] != -1)
358  {
359  verts.append(edgeToDualPoint_[edgeI]);
360  }
361  if (faceToDualPoint_[ie.faceLabel()] != -1)
362  {
363  doneEFaces[findIndex(eFaces, ie.faceLabel())] = true;
364  verts.append(faceToDualPoint_[ie.faceLabel()]);
365  }
366  if (cellToDualPoint_[ie.cellLabel()] != -1)
367  {
368  verts.append(cellToDualPoint_[ie.cellLabel()]);
369  }
370 
371  label currentDualCell0 = findDualCell(ie.cellLabel(), e[0]);
372  label currentDualCell1 = findDualCell(ie.cellLabel(), e[1]);
373 
374  ++ie;
375 
376  while (true)
377  {
378  label facei = ie.faceLabel();
379 
380  // Mark face as visited.
381  doneEFaces[findIndex(eFaces, facei)] = true;
382 
383  if (faceToDualPoint_[facei] != -1)
384  {
385  verts.append(faceToDualPoint_[facei]);
386  }
387 
388  label celli = ie.cellLabel();
389 
390  if (celli == -1)
391  {
392  // At ending boundary face. We've stored the face point above
393  // so this is the whole face.
394  break;
395  }
396 
397 
398  label dualCell0 = findDualCell(celli, e[0]);
399  label dualCell1 = findDualCell(celli, e[1]);
400 
401  // Generate face. (always if splitFace=true; only if needed to
402  // separate cells otherwise)
403  if
404  (
405  splitFace
406  || (
407  dualCell0 != currentDualCell0
408  || dualCell1 != currentDualCell1
409  )
410  )
411  {
412  // Close current face.
413  addInternalFace
414  (
415  -1, // masterPointi
416  edgeI, // masterEdgeI
417  -1, // masterFacei
418  edgeOrder,
419  currentDualCell0,
420  currentDualCell1,
421  verts.shrink(),
422  meshMod
423  );
424 
425  // Restart
426  currentDualCell0 = dualCell0;
427  currentDualCell1 = dualCell1;
428 
429  verts.clear();
430  if (edgeToDualPoint_[edgeI] != -1)
431  {
432  verts.append(edgeToDualPoint_[edgeI]);
433  }
434  if (faceToDualPoint_[facei] != -1)
435  {
436  verts.append(faceToDualPoint_[facei]);
437  }
438  }
439 
440  if (cellToDualPoint_[celli] != -1)
441  {
442  verts.append(cellToDualPoint_[celli]);
443  }
444 
445  ++ie;
446 
447  if (ie == ie.end())
448  {
449  // Back at start face (for internal edge only). See if this needs
450  // adding.
451  if (isBoundaryEdge.get(edgeI) == 0)
452  {
453  label startDual = faceToDualPoint_[startFaceLabel];
454 
455  if (startDual != -1 && findIndex(verts, startDual) == -1)
456  {
457  verts.append(startDual);
458  }
459  }
460  break;
461  }
462  }
463 
464  verts.shrink();
465  addInternalFace
466  (
467  -1, // masterPointi
468  edgeI, // masterEdgeI
469  -1, // masterFacei
470  edgeOrder,
471  currentDualCell0,
472  currentDualCell1,
473  verts,
474  meshMod
475  );
476 }
477 
478 
479 // Walks around circumference of facei. Creates single face. Gets given
480 // starting (feature) edge to start from. Returns ending edge. (all edges
481 // in form of index in faceEdges)
482 void Foam::meshDualiser::createFaceFromInternalFace
483 (
484  const label facei,
485  label& fp,
486  polyTopoChange& meshMod
487 ) const
488 {
489  const face& f = mesh_.faces()[facei];
490  const labelList& fEdges = mesh_.faceEdges()[facei];
491  label own = mesh_.faceOwner()[facei];
492  label nei = mesh_.faceNeighbour()[facei];
493 
494  // Pout<< "createFaceFromInternalFace : At face:" << facei
495  // << " verts:" << f
496  // << " points:" << UIndirectList<point>(mesh_.points(), f)()
497  // << " started walking at edge:" << fEdges[fp]
498  // << " verts:" << mesh_.edges()[fEdges[fp]]
499  // << endl;
500 
501 
502  // Walk and collect face.
503  DynamicList<label> verts(100);
504 
505  verts.append(faceToDualPoint_[facei]);
506  verts.append(edgeToDualPoint_[fEdges[fp]]);
507 
508  // Step to vertex after edge mid
509  fp = f.fcIndex(fp);
510 
511  // Get cells on either side of face at that point
512  label currentDualCell0 = findDualCell(own, f[fp]);
513  label currentDualCell1 = findDualCell(nei, f[fp]);
514 
515  forAll(f, i)
516  {
517  // Check vertex
518  if (pointToDualPoint_[f[fp]] != -1)
519  {
520  verts.append(pointToDualPoint_[f[fp]]);
521  }
522 
523  // Edge between fp and fp+1
524  label edgeI = fEdges[fp];
525 
526  if (edgeToDualPoint_[edgeI] != -1)
527  {
528  verts.append(edgeToDualPoint_[edgeI]);
529  }
530 
531  // Next vertex on edge
532  label nextFp = f.fcIndex(fp);
533 
534  // Get dual cells on nextFp to check whether face needs closing.
535  label dualCell0 = findDualCell(own, f[nextFp]);
536  label dualCell1 = findDualCell(nei, f[nextFp]);
537 
538  if (dualCell0 != currentDualCell0 || dualCell1 != currentDualCell1)
539  {
540  // Check: make sure that there is a midpoint on the edge.
541  if (edgeToDualPoint_[edgeI] == -1)
542  {
544  << "face:" << facei << " verts:" << f
545  << " points:" << UIndirectList<point>(mesh_.points(), f)()
546  << " no feature edge between " << f[fp]
547  << " and " << f[nextFp] << " although have different"
548  << " dual cells." << endl
549  << "point " << f[fp] << " has dual cells "
550  << currentDualCell0 << " and " << currentDualCell1
551  << " ; point "<< f[nextFp] << " has dual cells "
552  << dualCell0 << " and " << dualCell1
553  << abort(FatalError);
554  }
555 
556 
557  // Close current face.
558  verts.shrink();
559  addInternalFace
560  (
561  -1, // masterPointi
562  -1, // masterEdgeI
563  facei, // masterFacei
564  true, // edgeOrder,
565  currentDualCell0,
566  currentDualCell1,
567  verts,
568  meshMod
569  );
570  break;
571  }
572 
573  fp = nextFp;
574  }
575 }
576 
577 
578 // Given a point on a face converts the faces around the point.
579 // (pointFaces()). Gets starting face and marks off visited faces in donePFaces.
580 void Foam::meshDualiser::createFacesAroundBoundaryPoint
581 (
582  const label patchi,
583  const label patchPointi,
584  const label startFacei,
585  polyTopoChange& meshMod,
586  boolList& donePFaces // pFaces visited
587 ) const
588 {
589  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
590  const polyPatch& pp = patches[patchi];
591  const labelList& pFaces = pp.pointFaces()[patchPointi];
592  const labelList& own = mesh_.faceOwner();
593 
594  label pointi = pp.meshPoints()[patchPointi];
595 
596  if (pointToDualPoint_[pointi] == -1)
597  {
598  // Not a feature point. Loop over all connected
599  // pointFaces.
600 
601  // Starting face
602  label facei = startFacei;
603 
604  DynamicList<label> verts(4);
605 
606  while (true)
607  {
608  label index = findIndex(pFaces, facei-pp.start());
609 
610  // Has face been visited already?
611  if (donePFaces[index])
612  {
613  break;
614  }
615  donePFaces[index] = true;
616 
617  // Insert face centre
618  verts.append(faceToDualPoint_[facei]);
619 
620  label dualCelli = findDualCell(own[facei], pointi);
621 
622  // Get the edge before the patchPointi
623  const face& f = mesh_.faces()[facei];
624  label fp = findIndex(f, pointi);
625  label prevFp = f.rcIndex(fp);
626  label edgeI = mesh_.faceEdges()[facei][prevFp];
627 
628  if (edgeToDualPoint_[edgeI] != -1)
629  {
630  verts.append(edgeToDualPoint_[edgeI]);
631  }
632 
633  // Get next boundary face (whilst staying on edge).
634  edgeFaceCirculator circ
635  (
636  mesh_,
637  facei,
638  true, // ownerSide
639  prevFp, // index of edge in face
640  true // isBoundaryEdge
641  );
642 
643  do
644  {
645  ++circ;
646  }
647  while (mesh_.isInternalFace(circ.faceLabel()));
648 
649  // Step to next face
650  facei = circ.faceLabel();
651 
652  if (facei < pp.start() || facei >= pp.start()+pp.size())
653  {
655  << "Walked from face on patch:" << patchi
656  << " to face:" << facei
657  << " fc:" << mesh_.faceCentres()[facei]
658  << " on patch:" << patches.whichPatch(facei)
659  << abort(FatalError);
660  }
661 
662  // Check if different cell.
663  if (dualCelli != findDualCell(own[facei], pointi))
664  {
666  << "Different dual cells but no feature edge"
667  << " in between point:" << pointi
668  << " coord:" << mesh_.points()[pointi]
669  << abort(FatalError);
670  }
671  }
672 
673  verts.shrink();
674 
675  label dualCelli = findDualCell(own[facei], pointi);
676 
677  // Bit dodgy: create dualface from the last face (instead of from
678  // the central point). This will also use the original faceZone to
679  // put the new face (which might span multiple original faces) in.
680 
681  addBoundaryFace
682  (
683  // pointi, // masterPointi
684  -1, // masterPointi
685  -1, // masterEdgeI
686  facei, // masterFacei
687  dualCelli,
688  patchi,
689  verts,
690  meshMod
691  );
692  }
693  else
694  {
695  label facei = startFacei;
696 
697  // Storage for face
698  DynamicList<label> verts(mesh_.faces()[facei].size());
699 
700  // Starting point.
701  verts.append(pointToDualPoint_[pointi]);
702 
703  // Find edge between pointi and next point on face.
704  const labelList& fEdges = mesh_.faceEdges()[facei];
705  label nextEdgeI = fEdges[findIndex(mesh_.faces()[facei], pointi)];
706  if (edgeToDualPoint_[nextEdgeI] != -1)
707  {
708  verts.append(edgeToDualPoint_[nextEdgeI]);
709  }
710 
711  do
712  {
713  label index = findIndex(pFaces, facei-pp.start());
714 
715  // Has face been visited already?
716  if (donePFaces[index])
717  {
718  break;
719  }
720  donePFaces[index] = true;
721 
722  // Face centre
723  verts.append(faceToDualPoint_[facei]);
724 
725  // Find edge before pointi on facei
726  const labelList& fEdges = mesh_.faceEdges()[facei];
727  const face& f = mesh_.faces()[facei];
728  label prevFp = f.rcIndex(findIndex(f, pointi));
729  label edgeI = fEdges[prevFp];
730 
731  if (edgeToDualPoint_[edgeI] != -1)
732  {
733  // Feature edge. Close any face so far. Note: uses face to
734  // create dualFace from. Could use pointi instead.
735  verts.append(edgeToDualPoint_[edgeI]);
736  addBoundaryFace
737  (
738  -1, // masterPointi
739  -1, // masterEdgeI
740  facei, // masterFacei
741  findDualCell(own[facei], pointi),
742  patchi,
743  verts.shrink(),
744  meshMod
745  );
746  verts.clear();
747 
748  verts.append(pointToDualPoint_[pointi]);
749  verts.append(edgeToDualPoint_[edgeI]);
750  }
751 
752  // Cross edgeI to next boundary face
753  edgeFaceCirculator circ
754  (
755  mesh_,
756  facei,
757  true, // ownerSide
758  prevFp, // index of edge in face
759  true // isBoundaryEdge
760  );
761 
762  do
763  {
764  ++circ;
765  }
766  while (mesh_.isInternalFace(circ.faceLabel()));
767 
768  // Step to next face. Quit if not on same patch.
769  facei = circ.faceLabel();
770  }
771  while
772  (
773  facei != startFacei
774  && facei >= pp.start()
775  && facei < pp.start()+pp.size()
776  );
777 
778  if (verts.size() > 2)
779  {
780  // Note: face created from face, not from pointi
781  addBoundaryFace
782  (
783  -1, // masterPointi
784  -1, // masterEdgeI
785  startFacei, // masterFacei
786  findDualCell(own[facei], pointi),
787  patchi,
788  verts.shrink(),
789  meshMod
790  );
791  }
792  }
793 }
794 
795 
796 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
797 
798 Foam::meshDualiser::meshDualiser(const polyMesh& mesh)
799 :
800  mesh_(mesh),
801  pointToDualCells_(mesh_.nPoints()),
802  pointToDualPoint_(mesh_.nPoints(), -1),
803  cellToDualPoint_(mesh_.nCells()),
804  faceToDualPoint_(mesh_.nFaces(), -1),
805  edgeToDualPoint_(mesh_.nEdges(), -1)
806 {}
807 
808 
809 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
810 
812 (
813  const bool splitFace,
814  const labelList& featureFaces,
815  const labelList& featureEdges,
816  const labelList& singleCellFeaturePoints,
817  const labelList& multiCellFeaturePoints,
818  polyTopoChange& meshMod
819 )
820 {
821  const labelList& own = mesh_.faceOwner();
822  const labelList& nei = mesh_.faceNeighbour();
823  const vectorField& cellCentres = mesh_.cellCentres();
824 
825  // Mark boundary edges and points.
826  // (Note: in 1.4.2 we can use the built-in mesh point ordering
827  // facility instead)
828  PackedBoolList isBoundaryEdge(mesh_.nEdges());
829  for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
830  {
831  const labelList& fEdges = mesh_.faceEdges()[facei];
832 
833  forAll(fEdges, i)
834  {
835  isBoundaryEdge.set(fEdges[i], 1);
836  }
837  }
838 
839 
840  if (splitFace)
841  {
842  // This is a special mode where whenever we are walking around an edge
843  // every area through a cell becomes a separate dualface. So two
844  // dual cells will probably have more than one dualface between them!
845  // This mode implies that
846  // - all faces have to be feature faces since there has to be a
847  // dualpoint at the face centre.
848  // - all edges have to be feature edges ,,
849  boolList featureFaceSet(mesh_.nFaces(), false);
850  forAll(featureFaces, i)
851  {
852  featureFaceSet[featureFaces[i]] = true;
853  }
854  label facei = findIndex(featureFaceSet, false);
855 
856  if (facei != -1)
857  {
859  << "In split-face-mode (splitFace=true) but not all faces"
860  << " marked as feature faces." << endl
861  << "First conflicting face:" << facei
862  << " centre:" << mesh_.faceCentres()[facei]
863  << abort(FatalError);
864  }
865 
866  boolList featureEdgeSet(mesh_.nEdges(), false);
867  forAll(featureEdges, i)
868  {
869  featureEdgeSet[featureEdges[i]] = true;
870  }
871  label edgeI = findIndex(featureEdgeSet, false);
872 
873  if (edgeI != -1)
874  {
875  const edge& e = mesh_.edges()[edgeI];
877  << "In split-face-mode (splitFace=true) but not all edges"
878  << " marked as feature edges." << endl
879  << "First conflicting edge:" << edgeI
880  << " verts:" << e
881  << " coords:" << mesh_.points()[e[0]] << mesh_.points()[e[1]]
882  << abort(FatalError);
883  }
884  }
885  else
886  {
887  // Check that all boundary faces are feature faces.
888 
889  boolList featureFaceSet(mesh_.nFaces(), false);
890  forAll(featureFaces, i)
891  {
892  featureFaceSet[featureFaces[i]] = true;
893  }
894  for
895  (
896  label facei = mesh_.nInternalFaces();
897  facei < mesh_.nFaces();
898  facei++
899  )
900  {
901  if (!featureFaceSet[facei])
902  {
904  << "Not all boundary faces marked as feature faces."
905  << endl
906  << "First conflicting face:" << facei
907  << " centre:" << mesh_.faceCentres()[facei]
908  << abort(FatalError);
909  }
910  }
911  }
912 
913 
914 
915 
916  // Start creating cells, points, and faces (in that order)
917 
918 
919  // 1. Mark which cells to create
920  // Mostly every point becomes one cell but sometimes (for feature points)
921  // all cells surrounding a feature point become cells. Also a non-manifold
922  // point can create two cells! So a dual cell is uniquely defined by a
923  // mesh point + cell (as in pointCells index)
924 
925  // 2. Mark which face centres to create
926 
927  // 3. Internal faces can now consist of
928  // - only cell centres of walk around edge
929  // - cell centres + face centres of walk around edge
930  // - same but now other side is not a single cell
931 
932  // 4. Boundary faces (or internal faces between cell zones!) now consist of
933  // - walk around boundary point.
934 
935 
936 
937  autoPtr<OFstream> dualCcStr;
938  if (debug)
939  {
940  dualCcStr.reset(new OFstream("dualCc.obj"));
941  Pout<< "Dumping centres of dual cells to " << dualCcStr().name()
942  << endl;
943  }
944 
945 
946  // Dual cells (from points)
947  // ~~~~~~~~~~~~~~~~~~~~~~~~
948 
949  // pointToDualCells_[pointi]
950  // - single entry : all cells surrounding point all become the same
951  // cell.
952  // - multiple entries: in order of pointCells.
953 
954 
955  // feature points that become single cell
956  forAll(singleCellFeaturePoints, i)
957  {
958  label pointi = singleCellFeaturePoints[i];
959 
960  pointToDualPoint_[pointi] = meshMod.addPoint
961  (
962  mesh_.points()[pointi],
963  pointi, // masterPoint
964  mesh_.pointZones().whichZone(pointi), // zoneID
965  true // inCell
966  );
967 
968  // Generate single cell
969  pointToDualCells_[pointi].setSize(1);
970  pointToDualCells_[pointi][0] = meshMod.addCell
971  (
972  pointi, // masterPointID,
973  -1, // masterEdgeID,
974  -1, // masterFaceID,
975  -1, // masterCellID,
976  -1 // zoneID
977  );
978  if (dualCcStr.valid())
979  {
980  meshTools::writeOBJ(dualCcStr(), mesh_.points()[pointi]);
981  }
982  }
983 
984  // feature points that become multiple cells
985  forAll(multiCellFeaturePoints, i)
986  {
987  label pointi = multiCellFeaturePoints[i];
988 
989  if (pointToDualCells_[pointi].size() > 0)
990  {
992  << "Point " << pointi << " at:" << mesh_.points()[pointi]
993  << " is both in singleCellFeaturePoints"
994  << " and multiCellFeaturePoints."
995  << abort(FatalError);
996  }
997 
998  pointToDualPoint_[pointi] = meshMod.addPoint
999  (
1000  mesh_.points()[pointi],
1001  pointi, // masterPoint
1002  mesh_.pointZones().whichZone(pointi), // zoneID
1003  true // inCell
1004  );
1005 
1006  // Create dualcell for every cell connected to dual point
1007 
1008  const labelList& pCells = mesh_.pointCells()[pointi];
1009 
1010  pointToDualCells_[pointi].setSize(pCells.size());
1011 
1012  forAll(pCells, pCelli)
1013  {
1014  pointToDualCells_[pointi][pCelli] = meshMod.addCell
1015  (
1016  pointi, // masterPointID
1017  -1, // masterEdgeID
1018  -1, // masterFaceID
1019  -1, // masterCellID
1020  mesh_.cellZones().whichZone(pCells[pCelli]) // zoneID
1021  );
1022  if (dualCcStr.valid())
1023  {
1025  (
1026  dualCcStr(),
1027  0.5*(mesh_.points()[pointi]+cellCentres[pCells[pCelli]])
1028  );
1029  }
1030  }
1031  }
1032  // Normal points
1033  forAll(mesh_.points(), pointi)
1034  {
1035  if (pointToDualCells_[pointi].empty())
1036  {
1037  pointToDualCells_[pointi].setSize(1);
1038  pointToDualCells_[pointi][0] = meshMod.addCell
1039  (
1040  pointi, // masterPointID,
1041  -1, // masterEdgeID,
1042  -1, // masterFaceID,
1043  -1, // masterCellID,
1044  -1 // zoneID
1045  );
1046 
1047  if (dualCcStr.valid())
1048  {
1049  meshTools::writeOBJ(dualCcStr(), mesh_.points()[pointi]);
1050  }
1051  }
1052  }
1053 
1054 
1055  // Dual points (from cell centres, feature faces, feature edges)
1056  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1057 
1058  forAll(cellToDualPoint_, celli)
1059  {
1060  cellToDualPoint_[celli] = meshMod.addPoint
1061  (
1062  cellCentres[celli],
1063  mesh_.faces()[mesh_.cells()[celli][0]][0], // masterPoint
1064  -1, // zoneID
1065  true // inCell
1066  );
1067  }
1068 
1069  // From face to dual point
1070 
1071  forAll(featureFaces, i)
1072  {
1073  label facei = featureFaces[i];
1074 
1075  faceToDualPoint_[facei] = meshMod.addPoint
1076  (
1077  mesh_.faceCentres()[facei],
1078  mesh_.faces()[facei][0], // masterPoint
1079  -1, // zoneID
1080  true // inCell
1081  );
1082  }
1083  // Detect whether different dual cells on either side of a face. This
1084  // would necessitate having a dual face built from the face and thus a
1085  // dual point at the face centre.
1086  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
1087  {
1088  if (faceToDualPoint_[facei] == -1)
1089  {
1090  const face& f = mesh_.faces()[facei];
1091 
1092  forAll(f, fp)
1093  {
1094  label ownDualCell = findDualCell(own[facei], f[fp]);
1095  label neiDualCell = findDualCell(nei[facei], f[fp]);
1096 
1097  if (ownDualCell != neiDualCell)
1098  {
1099  faceToDualPoint_[facei] = meshMod.addPoint
1100  (
1101  mesh_.faceCentres()[facei],
1102  f[fp], // masterPoint
1103  -1, // zoneID
1104  true // inCell
1105  );
1106 
1107  break;
1108  }
1109  }
1110  }
1111  }
1112 
1113  // From edge to dual point
1114 
1115  forAll(featureEdges, i)
1116  {
1117  label edgeI = featureEdges[i];
1118 
1119  const edge& e = mesh_.edges()[edgeI];
1120 
1121  edgeToDualPoint_[edgeI] = meshMod.addPoint
1122  (
1123  e.centre(mesh_.points()),
1124  e[0], // masterPoint
1125  -1, // zoneID
1126  true // inCell
1127  );
1128  }
1129 
1130  // Detect whether different dual cells on either side of an edge. This
1131  // would necessitate having a dual face built perpendicular to the edge
1132  // and thus a dual point at the mid of the edge.
1133  // Note: not really true - the face can be built without the edge centre!
1134  const labelListList& edgeCells = mesh_.edgeCells();
1135 
1136  forAll(edgeCells, edgeI)
1137  {
1138  if (edgeToDualPoint_[edgeI] == -1)
1139  {
1140  const edge& e = mesh_.edges()[edgeI];
1141 
1142  // We need a point on the edge if not all cells on both sides
1143  // are the same.
1144 
1145  const labelList& eCells = mesh_.edgeCells()[edgeI];
1146 
1147  label dualE0 = findDualCell(eCells[0], e[0]);
1148  label dualE1 = findDualCell(eCells[0], e[1]);
1149 
1150  for (label i = 1; i < eCells.size(); i++)
1151  {
1152  label newDualE0 = findDualCell(eCells[i], e[0]);
1153 
1154  if (dualE0 != newDualE0)
1155  {
1156  edgeToDualPoint_[edgeI] = meshMod.addPoint
1157  (
1158  e.centre(mesh_.points()),
1159  e[0], // masterPoint
1160  mesh_.pointZones().whichZone(e[0]), // zoneID
1161  true // inCell
1162  );
1163 
1164  break;
1165  }
1166 
1167  label newDualE1 = findDualCell(eCells[i], e[1]);
1168 
1169  if (dualE1 != newDualE1)
1170  {
1171  edgeToDualPoint_[edgeI] = meshMod.addPoint
1172  (
1173  e.centre(mesh_.points()),
1174  e[1], // masterPoint
1175  mesh_.pointZones().whichZone(e[1]), // zoneID
1176  true // inCell
1177  );
1178 
1179  break;
1180  }
1181  }
1182  }
1183  }
1184 
1185  // Make sure all boundary edges emanating from feature points are
1186  // feature edges as well.
1187  forAll(singleCellFeaturePoints, i)
1188  {
1189  generateDualBoundaryEdges
1190  (
1191  isBoundaryEdge,
1192  singleCellFeaturePoints[i],
1193  meshMod
1194  );
1195  }
1196  forAll(multiCellFeaturePoints, i)
1197  {
1198  generateDualBoundaryEdges
1199  (
1200  isBoundaryEdge,
1201  multiCellFeaturePoints[i],
1202  meshMod
1203  );
1204  }
1205 
1206 
1207  // Check for duplicate points
1208  if (debug)
1209  {
1210  dumpPolyTopoChange(meshMod, "generatedPoints_");
1211  }
1212 
1213 
1214  // Now we have all points and cells
1215  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1216  // - pointToDualCells_ : per point a single dualCell or multiple dualCells
1217  // - pointToDualPoint_ : per point -1 or the dual point at the coordinate
1218  // - edgeToDualPoint_ : per edge -1 or the edge centre
1219  // - faceToDualPoint_ : per face -1 or the face centre
1220  // - cellToDualPoint_ : per cell the cell centre
1221  // Now we have to walk all edges and construct faces. Either single face
1222  // per edge or multiple (-if nonmanifold edge -if different dualcells)
1223 
1224  const edgeList& edges = mesh_.edges();
1225 
1226  forAll(edges, edgeI)
1227  {
1228  const labelList& eFaces = mesh_.edgeFaces()[edgeI];
1229 
1230  boolList doneEFaces(eFaces.size(), false);
1231 
1232  forAll(eFaces, i)
1233  {
1234  if (!doneEFaces[i])
1235  {
1236  // We found a face that hasn't yet been visited. This might
1237  // happen for non-manifold edges where a single edge can
1238  // become multiple faces.
1239 
1240  label startFacei = eFaces[i];
1241 
1242  // Pout<< "Walking edge:" << edgeI
1243  // << " points:" << mesh_.points()[e[0]]
1244  // << mesh_.points()[e[1]]
1245  // << " startFace:" << startFacei
1246  // << " at:" << mesh_.faceCentres()[startFacei]
1247  // << endl;
1248 
1249  createFacesAroundEdge
1250  (
1251  splitFace,
1252  isBoundaryEdge,
1253  edgeI,
1254  startFacei,
1255  meshMod,
1256  doneEFaces
1257  );
1258  }
1259  }
1260  }
1261 
1262  if (debug)
1263  {
1264  dumpPolyTopoChange(meshMod, "generatedFacesFromEdges_");
1265  }
1266 
1267  // Create faces from feature faces. These can be internal or external faces.
1268  // - feature face : centre needs to be included.
1269  // - single cells on either side: triangulate
1270  // - multiple cells: create single face between unique cell pair. Only
1271  // create face where cells differ on either side.
1272  // - non-feature face : in between cell zones.
1273  forAll(faceToDualPoint_, facei)
1274  {
1275  if (faceToDualPoint_[facei] != -1 && mesh_.isInternalFace(facei))
1276  {
1277  const face& f = mesh_.faces()[facei];
1278  const labelList& fEdges = mesh_.faceEdges()[facei];
1279 
1280  // Starting edge
1281  label fp = 0;
1282 
1283  do
1284  {
1285  // Find edge that is in dual mesh and where the
1286  // next point (fp+1) has different dual cells on either side.
1287  bool foundStart = false;
1288 
1289  do
1290  {
1291  if
1292  (
1293  edgeToDualPoint_[fEdges[fp]] != -1
1294  && !sameDualCell(facei, f.nextLabel(fp))
1295  )
1296  {
1297  foundStart = true;
1298  break;
1299  }
1300  fp = f.fcIndex(fp);
1301  }
1302  while (fp != 0);
1303 
1304  if (!foundStart)
1305  {
1306  break;
1307  }
1308 
1309  // Walk from edge fp and generate a face.
1310  createFaceFromInternalFace
1311  (
1312  facei,
1313  fp,
1314  meshMod
1315  );
1316  }
1317  while (fp != 0);
1318  }
1319  }
1320 
1321  if (debug)
1322  {
1323  dumpPolyTopoChange(meshMod, "generatedFacesFromFeatFaces_");
1324  }
1325 
1326 
1327  // Create boundary faces. Every boundary point has one or more dualcells.
1328  // These need to be closed.
1329  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1330 
1332  {
1333  const polyPatch& pp = patches[patchi];
1334 
1335  const labelListList& pointFaces = pp.pointFaces();
1336 
1337  forAll(pointFaces, patchPointi)
1338  {
1339  const labelList& pFaces = pointFaces[patchPointi];
1340 
1341  boolList donePFaces(pFaces.size(), false);
1342 
1343  forAll(pFaces, i)
1344  {
1345  if (!donePFaces[i])
1346  {
1347  // Starting face
1348  label startFacei = pp.start()+pFaces[i];
1349 
1350  // Pout<< "Walking around point:" << pointi
1351  // << " coord:" << mesh_.points()[pointi]
1352  // << " on patch:" << patchi
1353  // << " startFace:" << startFacei
1354  // << " at:" << mesh_.faceCentres()[startFacei]
1355  // << endl;
1356 
1357  createFacesAroundBoundaryPoint
1358  (
1359  patchi,
1360  patchPointi,
1361  startFacei,
1362  meshMod,
1363  donePFaces // pFaces visited
1364  );
1365  }
1366  }
1367  }
1368  }
1369 
1370  if (debug)
1371  {
1372  dumpPolyTopoChange(meshMod, "generatedFacesFromBndFaces_");
1373  }
1374 }
1375 
1376 
1377 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void setSize(const label)
Reset size of List.
Definition: List.C:281
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: MeshZones.C:221
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:85
label rcIndex(const label i) const
Return the reverse circular index, i.e. the previous index.
Definition: UListI.H:65
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
static label getMinIndex(const face &f, const label v0, const label v1)
Helper: find index in face of edge or -1. Index is such that edge is.
void setRefinement(const bool splitFace, const labelList &featureFaces, const labelList &featureEdges, const labelList &singleCellFeaturePoints, const labelList &multiCellFeaturePoints, polyTopoChange &meshMod)
Insert all changes into meshMod to convert the polyMesh into.
meshDualiser(const polyMesh &)
Construct from mesh.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1373
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1386
const meshPointZones & pointZones() const
Return point zones.
Definition: polyMesh.H:441
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1392
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1360
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:453
label nEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const vectorField & faceCentres() const
const labelListList & pointCells() const
const vectorField & cellCentres() const
label nInternalFaces() const
label nFaces() const
const labelListList & edgeFaces() const
const labelListList & faceEdges() const
const labelListList & edgeCells() const
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
const cellList & cells() const
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
const pointField & points
label nPoints
const fvPatchList & patches
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:105
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
defineTypeNameAndDebug(combustionModel, 0)
Field< vector > vectorField
Specialisation of Field<T> for vector.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
List< edge > edgeList
Definition: edgeList.H:38
static const char nl
Definition: Ostream.H:260
labelList f(nPoints)
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< small) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &mergedCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:229