meshCutter.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 "meshCutter.H"
27 #include "polyMesh.H"
28 #include "polyTopoChange.H"
29 #include "cellCuts.H"
30 #include "mapPolyMesh.H"
31 #include "meshTools.H"
32 #include "polyModifyFace.H"
33 #include "polyAddPoint.H"
34 #include "polyAddFace.H"
35 #include "polyAddCell.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(meshCutter, 0);
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Static Functions * * * * * * * * * * * //
46 
47 bool Foam::meshCutter::uses(const labelList& elems1, const labelList& elems2)
48 {
49  forAll(elems1, elemI)
50  {
51  if (findIndex(elems2, elems1[elemI]) != -1)
52  {
53  return true;
54  }
55  }
56  return false;
57 }
58 
59 
60 bool Foam::meshCutter::isIn
61 (
62  const edge& twoCuts,
63  const labelList& cuts
64 )
65 {
66  label index = findIndex(cuts, twoCuts[0]);
67 
68  if (index == -1)
69  {
70  return false;
71  }
72 
73  return
74  (
75  cuts[cuts.fcIndex(index)] == twoCuts[1]
76  || cuts[cuts.rcIndex(index)] == twoCuts[1]
77  );
78 }
79 
80 
81 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
82 
83 Foam::label Foam::meshCutter::findCutCell
84 (
85  const cellCuts& cuts,
86  const labelList& cellLabels
87 ) const
88 {
89  forAll(cellLabels, labelI)
90  {
91  label celli = cellLabels[labelI];
92 
93  if (cuts.cellLoops()[celli].size())
94  {
95  return celli;
96  }
97  }
98  return -1;
99 }
100 
101 
102 Foam::label Foam::meshCutter::findInternalFacePoint
103 (
104  const labelList& pointLabels
105 ) const
106 {
107  forAll(pointLabels, labelI)
108  {
109  label pointi = pointLabels[labelI];
110 
111  const labelList& pFaces = mesh().pointFaces()[pointi];
112 
113  forAll(pFaces, pFacei)
114  {
115  label facei = pFaces[pFacei];
116 
117  if (mesh().isInternalFace(facei))
118  {
119  return pointi;
120  }
121  }
122  }
123 
124  if (pointLabels.empty())
125  {
127  << "Empty pointLabels" << abort(FatalError);
128  }
129 
130  return -1;
131 }
132 
133 
134 void Foam::meshCutter::faceCells
135 (
136  const cellCuts& cuts,
137  const label facei,
138  label& own,
139  label& nei
140 ) const
141 {
142  const labelListList& anchorPts = cuts.cellAnchorPoints();
143  const labelListList& cellLoops = cuts.cellLoops();
144 
145  const face& f = mesh().faces()[facei];
146 
147  own = mesh().faceOwner()[facei];
148 
149  if (cellLoops[own].size() && uses(f, anchorPts[own]))
150  {
151  own = addedCells_[own];
152  }
153 
154  nei = -1;
155 
156  if (mesh().isInternalFace(facei))
157  {
158  nei = mesh().faceNeighbour()[facei];
159 
160  if (cellLoops[nei].size() && uses(f, anchorPts[nei]))
161  {
162  nei = addedCells_[nei];
163  }
164  }
165 }
166 
167 
168 void Foam::meshCutter::getFaceInfo
169 (
170  const label facei,
171  label& patchID,
172  label& zoneID,
173  label& zoneFlip
174 ) const
175 {
176  patchID = -1;
177 
178  if (!mesh().isInternalFace(facei))
179  {
180  patchID = mesh().boundaryMesh().whichPatch(facei);
181  }
182 
183  zoneID = mesh().faceZones().whichZone(facei);
184 
185  zoneFlip = false;
186 
187  if (zoneID >= 0)
188  {
189  const faceZone& fZone = mesh().faceZones()[zoneID];
190 
191  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
192  }
193 }
194 
195 
196 void Foam::meshCutter::addFace
197 (
198  polyTopoChange& meshMod,
199  const label facei,
200  const face& newFace,
201  const label own,
202  const label nei
203 )
204 {
205  label patchID, zoneID, zoneFlip;
206 
207  getFaceInfo(facei, patchID, zoneID, zoneFlip);
208 
209  if ((nei == -1) || (own < nei))
210  {
211  // Ordering ok.
212  if (debug & 2)
213  {
214  Pout<< "Adding face " << newFace
215  << " with new owner:" << own
216  << " with new neighbour:" << nei
217  << " patchID:" << patchID
218  << " zoneID:" << zoneID
219  << " zoneFlip:" << zoneFlip
220  << endl;
221  }
222 
223  meshMod.setAction
224  (
225  polyAddFace
226  (
227  newFace, // face
228  own, // owner
229  nei, // neighbour
230  -1, // master point
231  -1, // master edge
232  facei, // master face for addition
233  false, // flux flip
234  patchID, // patch for face
235  zoneID, // zone for face
236  zoneFlip // face zone flip
237  )
238  );
239  }
240  else
241  {
242  // Reverse owner/neighbour
243  if (debug & 2)
244  {
245  Pout<< "Adding (reversed) face " << newFace.reverseFace()
246  << " with new owner:" << nei
247  << " with new neighbour:" << own
248  << " patchID:" << patchID
249  << " zoneID:" << zoneID
250  << " zoneFlip:" << zoneFlip
251  << endl;
252  }
253 
254  meshMod.setAction
255  (
256  polyAddFace
257  (
258  newFace.reverseFace(), // face
259  nei, // owner
260  own, // neighbour
261  -1, // master point
262  -1, // master edge
263  facei, // master face for addition
264  false, // flux flip
265  patchID, // patch for face
266  zoneID, // zone for face
267  zoneFlip // face zone flip
268  )
269  );
270  }
271 }
272 
273 
274 void Foam::meshCutter::modFace
275 (
276  polyTopoChange& meshMod,
277  const label facei,
278  const face& newFace,
279  const label own,
280  const label nei
281 )
282 {
283  label patchID, zoneID, zoneFlip;
284 
285  getFaceInfo(facei, patchID, zoneID, zoneFlip);
286 
287  if
288  (
289  (own != mesh().faceOwner()[facei])
290  || (
291  mesh().isInternalFace(facei)
292  && (nei != mesh().faceNeighbour()[facei])
293  )
294  || (newFace != mesh().faces()[facei])
295  )
296  {
297  if (debug & 2)
298  {
299  Pout<< "Modifying face " << facei
300  << " old vertices:" << mesh().faces()[facei]
301  << " new vertices:" << newFace
302  << " new owner:" << own
303  << " new neighbour:" << nei
304  << " new zoneID:" << zoneID
305  << " new zoneFlip:" << zoneFlip
306  << endl;
307  }
308 
309  if ((nei == -1) || (own < nei))
310  {
311  meshMod.setAction
312  (
313  polyModifyFace
314  (
315  newFace, // modified face
316  facei, // label of face being modified
317  own, // owner
318  nei, // neighbour
319  false, // face flip
320  patchID, // patch for face
321  false, // remove from zone
322  zoneID, // zone for face
323  zoneFlip // face flip in zone
324  )
325  );
326  }
327  else
328  {
329  meshMod.setAction
330  (
331  polyModifyFace
332  (
333  newFace.reverseFace(), // modified face
334  facei, // label of face being modified
335  nei, // owner
336  own, // neighbour
337  false, // face flip
338  patchID, // patch for face
339  false, // remove from zone
340  zoneID, // zone for face
341  zoneFlip // face flip in zone
342  )
343  );
344  }
345  }
346 }
347 
348 
349 void Foam::meshCutter::copyFace
350 (
351  const face& f,
352  const label startFp,
353  const label endFp,
354  face& newFace
355 ) const
356 {
357  label fp = startFp;
358 
359  label newFp = 0;
360 
361  while (fp != endFp)
362  {
363  newFace[newFp++] = f[fp];
364 
365  fp = (fp + 1) % f.size();
366  }
367  newFace[newFp] = f[fp];
368 }
369 
370 
371 void Foam::meshCutter::splitFace
372 (
373  const face& f,
374  const label v0,
375  const label v1,
376 
377  face& f0,
378  face& f1
379 ) const
380 {
381  // Check if we find any new vertex which is part of the splitEdge.
382  label startFp = findIndex(f, v0);
383 
384  if (startFp == -1)
385  {
387  << "Cannot find vertex (new numbering) " << v0
388  << " on face " << f
389  << abort(FatalError);
390  }
391 
392  label endFp = findIndex(f, v1);
393 
394  if (endFp == -1)
395  {
397  << "Cannot find vertex (new numbering) " << v1
398  << " on face " << f
399  << abort(FatalError);
400  }
401 
402 
403  f0.setSize((endFp + 1 + f.size() - startFp) % f.size());
404  f1.setSize(f.size() - f0.size() + 2);
405 
406  copyFace(f, startFp, endFp, f0);
407  copyFace(f, endFp, startFp, f1);
408 }
409 
410 
411 Foam::face Foam::meshCutter::addEdgeCutsToFace(const label facei) const
412 {
413  const face& f = mesh().faces()[facei];
414 
415  face newFace(2 * f.size());
416 
417  label newFp = 0;
418 
419  forAll(f, fp)
420  {
421  // Duplicate face vertex .
422  newFace[newFp++] = f[fp];
423 
424  // Check if edge has been cut.
425  label fp1 = f.fcIndex(fp);
426 
427  HashTable<label, edge, Hash<edge>>::const_iterator fnd =
428  addedPoints_.find(edge(f[fp], f[fp1]));
429 
430  if (fnd != addedPoints_.end())
431  {
432  // edge has been cut. Introduce new vertex.
433  newFace[newFp++] = fnd();
434  }
435  }
436 
437  newFace.setSize(newFp);
438 
439  return newFace;
440 }
441 
442 
443 Foam::face Foam::meshCutter::loopToFace
444 (
445  const label celli,
446  const labelList& loop
447 ) const
448 {
449  face newFace(2*loop.size());
450 
451  label newFacei = 0;
452 
453  forAll(loop, fp)
454  {
455  label cut = loop[fp];
456 
457  if (isEdge(cut))
458  {
459  label edgeI = getEdge(cut);
460 
461  const edge& e = mesh().edges()[edgeI];
462 
463  label vertI = addedPoints_[e];
464 
465  newFace[newFacei++] = vertI;
466  }
467  else
468  {
469  // cut is vertex.
470  label vertI = getVertex(cut);
471 
472  newFace[newFacei++] = vertI;
473 
474  label nextCut = loop[loop.fcIndex(fp)];
475 
476  if (!isEdge(nextCut))
477  {
478  // From vertex to vertex -> cross cut only if no existing edge.
479  label nextVertI = getVertex(nextCut);
480 
481  label edgeI = meshTools::findEdge(mesh(), vertI, nextVertI);
482 
483  if (edgeI != -1)
484  {
485  // Existing edge. Insert split-edge point if any.
486  HashTable<label, edge, Hash<edge>>::const_iterator fnd =
487  addedPoints_.find(mesh().edges()[edgeI]);
488 
489  if (fnd != addedPoints_.end())
490  {
491  newFace[newFacei++] = fnd();
492  }
493  }
494  }
495  }
496  }
497  newFace.setSize(newFacei);
498 
499  return newFace;
500 }
501 
502 
503 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
504 
505 Foam::meshCutter::meshCutter(const polyMesh& mesh)
506 :
507  edgeVertex(mesh),
508  addedCells_(),
509  addedFaces_(),
510  addedPoints_()
511 
512 {}
513 
514 
515 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
516 
518 {}
519 
520 
521 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
522 
524 (
525  const cellCuts& cuts,
526  polyTopoChange& meshMod
527 )
528 {
529  // Clear and size maps here since mesh size will change.
530  addedCells_.clear();
531  addedCells_.resize(cuts.nLoops());
532 
533  addedFaces_.clear();
534  addedFaces_.resize(cuts.nLoops());
535 
536  addedPoints_.clear();
537  addedPoints_.resize(cuts.nLoops());
538 
539  if (cuts.nLoops() == 0)
540  {
541  return;
542  }
543 
544  const labelListList& anchorPts = cuts.cellAnchorPoints();
545  const labelListList& cellLoops = cuts.cellLoops();
546 
547  //
548  // Add new points along cut edges.
549  //
550 
551  forAll(cuts.edgeIsCut(), edgeI)
552  {
553  if (cuts.edgeIsCut()[edgeI])
554  {
555  const edge& e = mesh().edges()[edgeI];
556 
557  // Check if there is any cell using this edge.
558  if (debug && findCutCell(cuts, mesh().edgeCells()[edgeI]) == -1)
559  {
561  << "Problem: cut edge but none of the cells using it is\n"
562  << "edge:" << edgeI << " verts:" << e
563  << abort(FatalError);
564  }
565 
566  // One of the edge end points should be master point of nbCelli.
567  label masterPointi = e.start();
568 
569  const point& v0 = mesh().points()[e.start()];
570  const point& v1 = mesh().points()[e.end()];
571 
572  scalar weight = cuts.edgeWeight()[edgeI];
573 
574  point newPt = weight*v1 + (1.0-weight)*v0;
575 
576  label addedPointi =
577  meshMod.setAction
578  (
580  (
581  newPt, // point
582  masterPointi, // master point
583  -1, // zone for point
584  true // supports a cell
585  )
586  );
587 
588  // Store on (hash of) edge.
589  addedPoints_.insert(e, addedPointi);
590 
591  if (debug & 2)
592  {
593  Pout<< "Added point " << addedPointi
594  << " to vertex "
595  << masterPointi << " of edge " << edgeI
596  << " vertices " << e << endl;
597  }
598  }
599  }
600 
601  //
602  // Add cells (on 'anchor' side of cell)
603  //
604 
605  forAll(cellLoops, celli)
606  {
607  if (cellLoops[celli].size())
608  {
609  // Add a cell to the existing cell
610  label addedCelli =
611  meshMod.setAction
612  (
614  (
615  -1, // master point
616  -1, // master edge
617  -1, // master face
618  celli, // master cell
619  mesh().cellZones().whichZone(celli) // zone for cell
620  )
621  );
622 
623  addedCells_.insert(celli, addedCelli);
624 
625  if (debug & 2)
626  {
627  Pout<< "Added cell " << addedCells_[celli] << " to cell "
628  << celli << endl;
629  }
630  }
631  }
632 
633 
634  //
635  // For all cut cells add an internal face
636  //
637 
638  forAll(cellLoops, celli)
639  {
640  const labelList& loop = cellLoops[celli];
641 
642  if (loop.size())
643  {
644  // Convert loop (=list of cuts) into proper face.
645  // Orientation should already be ok. (done by cellCuts)
646  //
647  face newFace(loopToFace(celli, loop));
648 
649  // Pick any anchor point on cell
650  label masterPointi = findInternalFacePoint(anchorPts[celli]);
651 
652  label addedFacei =
653  meshMod.setAction
654  (
656  (
657  newFace, // face
658  celli, // owner
659  addedCells_[celli], // neighbour
660  masterPointi, // master point
661  -1, // master edge
662  -1, // master face for addition
663  false, // flux flip
664  -1, // patch for face
665  -1, // zone for face
666  false // face zone flip
667  )
668  );
669 
670  addedFaces_.insert(celli, addedFacei);
671 
672  if (debug & 2)
673  {
674  // Gets edgeweights of loop
675  scalarField weights(loop.size());
676  forAll(loop, i)
677  {
678  label cut = loop[i];
679 
680  weights[i] =
681  (
682  isEdge(cut)
683  ? cuts.edgeWeight()[getEdge(cut)]
684  : -GREAT
685  );
686  }
687 
688  Pout<< "Added splitting face " << newFace << " index:"
689  << addedFacei
690  << " to owner " << celli
691  << " neighbour " << addedCells_[celli]
692  << " from Loop:";
693  writeCuts(Pout, loop, weights);
694  Pout<< endl;
695  }
696  }
697  }
698 
699 
700  //
701  // Modify faces on the outside and create new ones
702  // (in effect split old faces into two)
703  //
704 
705  // Maintain whether face has been updated (for -split edges
706  // -new owner/neighbour)
707  boolList faceUptodate(mesh().nFaces(), false);
708 
709  const Map<edge>& faceSplitCuts = cuts.faceSplitCut();
710 
711  forAllConstIter(Map<edge>, faceSplitCuts, iter)
712  {
713  label facei = iter.key();
714 
715  // Renumber face to include split edges.
716  face newFace(addEdgeCutsToFace(facei));
717 
718  // Edge splitting the face. Convert cuts to new vertex numbering.
719  const edge& splitEdge = iter();
720 
721  label cut0 = splitEdge[0];
722 
723  label v0;
724  if (isEdge(cut0))
725  {
726  label edgeI = getEdge(cut0);
727  v0 = addedPoints_[mesh().edges()[edgeI]];
728  }
729  else
730  {
731  v0 = getVertex(cut0);
732  }
733 
734  label cut1 = splitEdge[1];
735  label v1;
736  if (isEdge(cut1))
737  {
738  label edgeI = getEdge(cut1);
739  v1 = addedPoints_[mesh().edges()[edgeI]];
740  }
741  else
742  {
743  v1 = getVertex(cut1);
744  }
745 
746  // Split face along the elements of the splitEdge.
747  face f0, f1;
748  splitFace(newFace, v0, v1, f0, f1);
749 
750  label own = mesh().faceOwner()[facei];
751 
752  label nei = -1;
753 
754  if (mesh().isInternalFace(facei))
755  {
756  nei = mesh().faceNeighbour()[facei];
757  }
758 
759  if (debug & 2)
760  {
761  Pout<< "Split face " << mesh().faces()[facei]
762  << " own:" << own << " nei:" << nei
763  << " into f0:" << f0
764  << " and f1:" << f1 << endl;
765  }
766 
767  // Check which face uses anchorPoints (connects to addedCell)
768  // and which one doesn't (connects to original cell)
769 
770  // Bit tricky. We have to know whether this faceSplit splits owner/
771  // neighbour or both. Even if cell is cut we have to make sure this is
772  // the one that cuts it (this face cut might not be the one splitting
773  // the cell)
774 
775  const face& f = mesh().faces()[facei];
776 
777  label f0Owner = -1;
778  label f1Owner = -1;
779 
780  if (cellLoops[own].empty())
781  {
782  f0Owner = own;
783  f1Owner = own;
784  }
785  else if (isIn(splitEdge, cellLoops[own]))
786  {
787  // Owner is cut by this splitCut. See which of f0, f1 gets
788  // owner, which gets addedCells_[owner]
789  if (uses(f0, anchorPts[own]))
790  {
791  f0Owner = addedCells_[own];
792  f1Owner = own;
793  }
794  else
795  {
796  f0Owner = own;
797  f1Owner = addedCells_[own];
798  }
799  }
800  else
801  {
802  // Owner not cut by this splitCut. Check on original face whether
803  // use anchorPts.
804  if (uses(f, anchorPts[own]))
805  {
806  label newCelli = addedCells_[own];
807  f0Owner = newCelli;
808  f1Owner = newCelli;
809  }
810  else
811  {
812  f0Owner = own;
813  f1Owner = own;
814  }
815  }
816 
817 
818  label f0Neighbour = -1;
819  label f1Neighbour = -1;
820 
821  if (nei != -1)
822  {
823  if (cellLoops[nei].empty())
824  {
825  f0Neighbour = nei;
826  f1Neighbour = nei;
827  }
828  else if (isIn(splitEdge, cellLoops[nei]))
829  {
830  // Neighbour is cut by this splitCut. See which of f0, f1
831  // gets which neighbour/addedCells_[neighbour]
832  if (uses(f0, anchorPts[nei]))
833  {
834  f0Neighbour = addedCells_[nei];
835  f1Neighbour = nei;
836  }
837  else
838  {
839  f0Neighbour = nei;
840  f1Neighbour = addedCells_[nei];
841  }
842  }
843  else
844  {
845  // neighbour not cut by this splitCut. Check on original face
846  // whether use anchorPts.
847  if (uses(f, anchorPts[nei]))
848  {
849  label newCelli = addedCells_[nei];
850  f0Neighbour = newCelli;
851  f1Neighbour = newCelli;
852  }
853  else
854  {
855  f0Neighbour = nei;
856  f1Neighbour = nei;
857  }
858  }
859  }
860 
861  // f0 is the added face, f1 the modified one
862  addFace(meshMod, facei, f0, f0Owner, f0Neighbour);
863 
864  modFace(meshMod, facei, f1, f1Owner, f1Neighbour);
865 
866  faceUptodate[facei] = true;
867  }
868 
869 
870  //
871  // Faces that have not been split but just appended to. Are guaranteed
872  // to be reachable from an edgeCut.
873  //
874 
875  const boolList& edgeIsCut = cuts.edgeIsCut();
876 
877  forAll(edgeIsCut, edgeI)
878  {
879  if (edgeIsCut[edgeI])
880  {
881  const labelList& eFaces = mesh().edgeFaces()[edgeI];
882 
883  forAll(eFaces, i)
884  {
885  label facei = eFaces[i];
886 
887  if (!faceUptodate[facei])
888  {
889  // Renumber face to include split edges.
890  face newFace(addEdgeCutsToFace(facei));
891 
892  if (debug & 2)
893  {
894  Pout<< "Added edge cuts to face " << facei
895  << " f:" << mesh().faces()[facei]
896  << " newFace:" << newFace << endl;
897  }
898 
899  // Get (new or original) owner and neighbour of facei
900  label own, nei;
901  faceCells(cuts, facei, own, nei);
902 
903  modFace(meshMod, facei, newFace, own, nei);
904 
905  faceUptodate[facei] = true;
906  }
907  }
908  }
909  }
910 
911 
912  //
913  // Correct any original faces on split cell for new neighbour/owner
914  //
915 
916  forAll(cellLoops, celli)
917  {
918  if (cellLoops[celli].size())
919  {
920  const labelList& cllFaces = mesh().cells()[celli];
921 
922  forAll(cllFaces, cllFacei)
923  {
924  label facei = cllFaces[cllFacei];
925 
926  if (!faceUptodate[facei])
927  {
928  // Update face with new owner/neighbour (if any)
929  const face& f = mesh().faces()[facei];
930 
931  if (debug && (f != addEdgeCutsToFace(facei)))
932  {
934  << "Problem: edges added to face which does "
935  << " not use a marked cut" << endl
936  << "facei:" << facei << endl
937  << "face:" << f << endl
938  << "newFace:" << addEdgeCutsToFace(facei)
939  << abort(FatalError);
940  }
941 
942  // Get (new or original) owner and neighbour of facei
943  label own, nei;
944  faceCells(cuts, facei, own, nei);
945 
946  modFace
947  (
948  meshMod,
949  facei,
950  f,
951  own,
952  nei
953  );
954 
955  faceUptodate[facei] = true;
956  }
957  }
958  }
959  }
960 
961  if (debug)
962  {
963  Pout<< "meshCutter:" << nl
964  << " cells split:" << addedCells_.size() << nl
965  << " faces added:" << addedFaces_.size() << nl
966  << " points added on edges:" << addedPoints_.size() << nl
967  << endl;
968  }
969 }
970 
971 
973 {
974  // Update stored labels for mesh change.
975 
976  {
977  // Create copy since new label might (temporarily) clash with existing
978  // key.
979  Map<label> newAddedCells(addedCells_.size());
980 
981  forAllConstIter(Map<label>, addedCells_, iter)
982  {
983  label celli = iter.key();
984  label newCelli = morphMap.reverseCellMap()[celli];
985 
986  label addedCelli = iter();
987 
988  label newAddedCelli = morphMap.reverseCellMap()[addedCelli];
989 
990  if (newCelli >= 0 && newAddedCelli >= 0)
991  {
992  if
993  (
994  (debug & 2)
995  && (newCelli != celli || newAddedCelli != addedCelli)
996  )
997  {
998  Pout<< "meshCutter::updateMesh :"
999  << " updating addedCell for cell " << celli
1000  << " from " << addedCelli
1001  << " to " << newAddedCelli << endl;
1002  }
1003  newAddedCells.insert(newCelli, newAddedCelli);
1004  }
1005  }
1006 
1007  // Copy
1008  addedCells_.transfer(newAddedCells);
1009  }
1010 
1011  {
1012  Map<label> newAddedFaces(addedFaces_.size());
1013 
1014  forAllConstIter(Map<label>, addedFaces_, iter)
1015  {
1016  label celli = iter.key();
1017  label newCelli = morphMap.reverseCellMap()[celli];
1018 
1019  label addedFacei = iter();
1020 
1021  label newAddedFacei = morphMap.reverseFaceMap()[addedFacei];
1022 
1023  if ((newCelli >= 0) && (newAddedFacei >= 0))
1024  {
1025  if
1026  (
1027  (debug & 2)
1028  && (newCelli != celli || newAddedFacei != addedFacei)
1029  )
1030  {
1031  Pout<< "meshCutter::updateMesh :"
1032  << " updating addedFace for cell " << celli
1033  << " from " << addedFacei
1034  << " to " << newAddedFacei
1035  << endl;
1036  }
1037  newAddedFaces.insert(newCelli, newAddedFacei);
1038  }
1039  }
1040 
1041  // Copy
1042  addedFaces_.transfer(newAddedFaces);
1043  }
1044 
1045  {
1046  HashTable<label, edge, Hash<edge>> newAddedPoints(addedPoints_.size());
1047 
1048  for
1049  (
1050  HashTable<label, edge, Hash<edge>>::const_iterator iter =
1051  addedPoints_.begin();
1052  iter != addedPoints_.end();
1053  ++iter
1054  )
1055  {
1056  const edge& e = iter.key();
1057 
1058  label newStart = morphMap.reversePointMap()[e.start()];
1059 
1060  label newEnd = morphMap.reversePointMap()[e.end()];
1061 
1062  label addedPointi = iter();
1063 
1064  label newAddedPointi = morphMap.reversePointMap()[addedPointi];
1065 
1066  if ((newStart >= 0) && (newEnd >= 0) && (newAddedPointi >= 0))
1067  {
1068  edge newE = edge(newStart, newEnd);
1069 
1070  if
1071  (
1072  (debug & 2)
1073  && (e != newE || newAddedPointi != addedPointi)
1074  )
1075  {
1076  Pout<< "meshCutter::updateMesh :"
1077  << " updating addedPoints for edge " << e
1078  << " from " << addedPointi
1079  << " to " << newAddedPointi
1080  << endl;
1081  }
1082 
1083  newAddedPoints.insert(newE, newAddedPointi);
1084  }
1085  }
1086 
1087  // Copy
1088  addedPoints_.transfer(newAddedPoints);
1089  }
1090 }
1091 
1092 
1093 // ************************************************************************* //
const labelListList & pointFaces() const
label end() const
Return end vertex label.
Definition: edgeI.H:92
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
const boolList & edgeIsCut() const
Is edge cut.
Definition: cellCuts.H:544
#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
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const double e
Elementary charge.
Definition: doubleFloat.H:78
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
label findEdge(const edgeList &edges, const labelList &candidates, const label v0, const label v1)
Return edge among candidates that uses the two vertices.
Definition: meshTools.C:337
Class containing data for point addition.
Definition: polyAddPoint.H:47
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
void transfer(HashTable< T, Key, Hash > &)
Transfer the contents of the argument table into this table.
Definition: HashTable.C:509
label nLoops() const
Number of valid cell loops.
Definition: cellCuts.H:579
Description of cuts across cells.
Definition: cellCuts.H:108
const cellList & cells() const
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
scalar f1
Definition: createFields.H:28
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:51
const labelListList & edgeFaces() const
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: meshCutter.C:972
Class containing data for cell addition.
Definition: polyAddCell.H:46
const polyMesh & mesh() const
Definition: edgeVertex.H:98
static const labelSphericalTensor labelI(1)
Identity labelTensor.
dynamicFvMesh & mesh
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:495
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
const scalarField & edgeWeight() const
If edge is cut gives weight (ratio between start() and end())
Definition: cellCuts.H:550
void clear()
Clear all entries from table.
Definition: HashTable.C:464
Combines edge or vertex in single label. Used to specify cuts across cell circumference.
Definition: edgeVertex.H:52
static label getVertex(const primitiveMesh &mesh, const label eVert)
Convert eVert to vertex label.
Definition: edgeVertex.H:142
List< label > labelList
A List of labels.
Definition: labelList.H:56
An STL-conforming hash table.
Definition: HashTable.H:61
static label getEdge(const primitiveMesh &mesh, const label eVert)
Convert eVert to edge label.
Definition: edgeVertex.H:126
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void setRefinement(const cellCuts &cuts, polyTopoChange &meshMod)
Do actual cutting with cut description. Inserts mesh changes.
Definition: meshCutter.C:524
const labelListList & cellLoops() const
For each cut cell the cut along the circumference.
Definition: cellCuts.H:573
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
label start() const
Return start vertex label.
Definition: edgeI.H:81
static const char nl
Definition: Ostream.H:262
defineTypeNameAndDebug(combustionModel, 0)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
static bool isEdge(const primitiveMesh &mesh, const label eVert)
Is eVert an edge?
Definition: edgeVertex.H:107
const Map< edge > & faceSplitCut() const
Gives for split face the two cuts that split the face into two.
Definition: cellCuts.H:567
const labelListList & cellAnchorPoints() const
For each cut cell the points on the &#39;anchor&#39; side of the cell.
Definition: cellCuts.H:585
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
~meshCutter()
Destructor.
Definition: meshCutter.C:517
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:221
Direct mesh changes based on v1.3 polyTopoChange syntax.
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:526
Hash function class for primitives. All non-primitives used to hash entries on hash tables likely nee...
Definition: Hash.H:54
void resize(const label newSize)
Resize the hash table for efficiency.
Definition: HashTable.C:428
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
Ostream & writeCuts(Ostream &os, const labelList &, const scalarField &) const
Write cut descriptions to Ostream.
Definition: edgeVertex.C:239
Namespace for OpenFOAM.
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:463
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.