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