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