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