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