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