cellCuts.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 "cellCuts.H"
27 #include "polyMesh.H"
28 #include "Time.H"
29 #include "ListOps.H"
30 #include "cellLooper.H"
31 #include "refineCell.H"
32 #include "meshTools.H"
33 #include "geomCellLooper.H"
34 #include "OFstream.H"
35 #include "plane.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(cellCuts, 0);
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Static Functions * * * * * * * * * * * //
46 
47 Foam::label Foam::cellCuts::findPartIndex
48 (
49  const labelList& elems,
50  const label nElems,
51  const label val
52 )
53 {
54  for (label i = 0; i < nElems; i++)
55  {
56  if (elems[i] == val)
57  {
58  return i;
59  }
60  }
61  return -1;
62 }
63 
64 
65 Foam::boolList Foam::cellCuts::expand
66 (
67  const label size,
68  const labelList& labels
69 )
70 {
71  boolList result(size, false);
72 
73  forAll(labels, labelI)
74  {
75  result[labels[labelI]] = true;
76  }
77  return result;
78 }
79 
80 
81 Foam::scalarField Foam::cellCuts::expand
82 (
83  const label size,
84  const labelList& labels,
85  const scalarField& weights
86 )
87 {
88  scalarField result(size, -GREAT);
89 
90  forAll(labels, labelI)
91  {
92  result[labels[labelI]] = weights[labelI];
93  }
94  return result;
95 }
96 
97 
98 Foam::label Foam::cellCuts::firstUnique
99 (
100  const labelList& lst,
101  const Map<label>& map
102 )
103 {
104  forAll(lst, i)
105  {
106  if (!map.found(lst[i]))
107  {
108  return i;
109  }
110  }
111  return -1;
112 }
113 
114 
115 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
116 
117 void Foam::cellCuts::writeUncutOBJ
118 (
119  const fileName& dir,
120  const label celli
121 ) const
122 {
123  // Cell edges
124  OFstream cutsStream(dir / "cell_" + name(celli) + ".obj");
125 
126  Pout<< "Writing cell for time " << mesh().time().timeName()
127  << " to " << cutsStream.name() << nl;
128 
130  (
131  cutsStream,
132  mesh().cells(),
133  mesh().faces(),
134  mesh().points(),
135  labelList(1, celli)
136  );
137 
138  // Loop cutting cell in two
139  OFstream cutStream(dir / "cellCuts_" + name(celli) + ".obj");
140 
141  Pout<< "Writing raw cuts on cell for time " << mesh().time().timeName()
142  << " to " << cutStream.name() << nl;
143 
144  const labelList& cPoints = mesh().cellPoints()[celli];
145 
146  forAll(cPoints, i)
147  {
148  label pointi = cPoints[i];
149  if (pointIsCut_[pointi])
150  {
151  meshTools::writeOBJ(cutStream, mesh().points()[pointi]);
152  }
153  }
154 
155  const pointField& pts = mesh().points();
156 
157  const labelList& cEdges = mesh().cellEdges()[celli];
158 
159  forAll(cEdges, i)
160  {
161  label edgeI = cEdges[i];
162 
163  if (edgeIsCut_[edgeI])
164  {
165  const edge& e = mesh().edges()[edgeI];
166 
167  const scalar w = edgeWeight_[edgeI];
168 
169  meshTools::writeOBJ(cutStream, w*pts[e[1]] + (1-w)*pts[e[0]]);
170  }
171  }
172 }
173 
174 
175 void Foam::cellCuts::writeOBJ
176 (
177  const fileName& dir,
178  const label celli,
179  const pointField& loopPoints,
180  const labelList& anchors
181 ) const
182 {
183  // Cell edges
184  OFstream cutsStream(dir / "cell_" + name(celli) + ".obj");
185 
186  Pout<< "Writing cell for time " << mesh().time().timeName()
187  << " to " << cutsStream.name() << nl;
188 
190  (
191  cutsStream,
192  mesh().cells(),
193  mesh().faces(),
194  mesh().points(),
195  labelList(1, celli)
196  );
197 
198 
199  // Loop cutting cell in two
200  OFstream loopStream(dir / "cellLoop_" + name(celli) + ".obj");
201 
202  Pout<< "Writing loop for time " << mesh().time().timeName()
203  << " to " << loopStream.name() << nl;
204 
205  label vertI = 0;
206 
207  writeOBJ(loopStream, loopPoints, vertI);
208 
209 
210  // Anchors for cell
211  OFstream anchorStream(dir / "anchors_" + name(celli) + ".obj");
212 
213  Pout<< "Writing anchors for time " << mesh().time().timeName()
214  << " to " << anchorStream.name() << endl;
215 
216  forAll(anchors, i)
217  {
218  meshTools::writeOBJ(anchorStream, mesh().points()[anchors[i]]);
219  }
220 }
221 
222 
223 Foam::label Foam::cellCuts::edgeEdgeToFace
224 (
225  const label celli,
226  const label edgeA,
227  const label edgeB
228 ) const
229 {
230  const labelList& cFaces = mesh().cells()[celli];
231 
232  forAll(cFaces, cFacei)
233  {
234  label facei = cFaces[cFacei];
235 
236  const labelList& fEdges = mesh().faceEdges()[facei];
237 
238  if
239  (
240  findIndex(fEdges, edgeA) != -1
241  && findIndex(fEdges, edgeB) != -1
242  )
243  {
244  return facei;
245  }
246  }
247 
248  // Coming here means the loop is illegal since the two edges
249  // are not shared by a face. We just mark loop as invalid and continue.
250 
252  << "cellCuts : Cannot find face on cell "
253  << celli << " that has both edges " << edgeA << ' ' << edgeB << endl
254  << "faces : " << cFaces << endl
255  << "edgeA : " << mesh().edges()[edgeA] << endl
256  << "edgeB : " << mesh().edges()[edgeB] << endl
257  << "Marking the loop across this cell as invalid" << endl;
258 
259  return -1;
260 }
261 
262 
263 Foam::label Foam::cellCuts::edgeVertexToFace
264 (
265  const label celli,
266  const label edgeI,
267  const label vertI
268 ) const
269 {
270  const labelList& cFaces = mesh().cells()[celli];
271 
272  forAll(cFaces, cFacei)
273  {
274  label facei = cFaces[cFacei];
275 
276  const face& f = mesh().faces()[facei];
277 
278  const labelList& fEdges = mesh().faceEdges()[facei];
279 
280  if
281  (
282  findIndex(fEdges, edgeI) != -1
283  && findIndex(f, vertI) != -1
284  )
285  {
286  return facei;
287  }
288  }
289 
291  << "cellCuts : Cannot find face on cell "
292  << celli << " that has both edge " << edgeI << " and vertex "
293  << vertI << endl
294  << "faces : " << cFaces << endl
295  << "edge : " << mesh().edges()[edgeI] << endl
296  << "Marking the loop across this cell as invalid" << endl;
297 
298  return -1;
299 }
300 
301 
302 Foam::label Foam::cellCuts::vertexVertexToFace
303 (
304  const label celli,
305  const label vertA,
306  const label vertB
307 ) const
308 {
309  const labelList& cFaces = mesh().cells()[celli];
310 
311  forAll(cFaces, cFacei)
312  {
313  label facei = cFaces[cFacei];
314 
315  const face& f = mesh().faces()[facei];
316 
317  if
318  (
319  findIndex(f, vertA) != -1
320  && findIndex(f, vertB) != -1
321  )
322  {
323  return facei;
324  }
325  }
326 
328  << "cellCuts : Cannot find face on cell "
329  << celli << " that has vertex " << vertA << " and vertex "
330  << vertB << endl
331  << "faces : " << cFaces << endl
332  << "Marking the loop across this cell as invalid" << endl;
333 
334  return -1;
335 }
336 
337 
338 void Foam::cellCuts::calcFaceCuts() const
339 {
340  if (faceCutsPtr_)
341  {
343  << "faceCuts already calculated" << abort(FatalError);
344  }
345 
346  const faceList& faces = mesh().faces();
347 
348 
349  faceCutsPtr_ = new labelListList(mesh().nFaces());
350  labelListList& faceCuts = *faceCutsPtr_;
351 
352  for (label facei = 0; facei < mesh().nFaces(); facei++)
353  {
354  const face& f = faces[facei];
355 
356  // Big enough storage (possibly all points and all edges cut). Shrink
357  // later on.
358  labelList& cuts = faceCuts[facei];
359 
360  cuts.setSize(2*f.size());
361 
362  label cutI = 0;
363 
364  // Do point-edge-point walk over face and collect all cuts.
365  // Problem is that we want to start from one of the endpoints of a
366  // string of connected cuts; we don't want to start somewhere in the
367  // middle.
368 
369  // Pass1: find first point cut not preceeded by a cut.
370  label startFp = -1;
371 
372  forAll(f, fp)
373  {
374  label v0 = f[fp];
375 
376  if (pointIsCut_[v0])
377  {
378  label vMin1 = f[f.rcIndex(fp)];
379 
380  if
381  (
382  !pointIsCut_[vMin1]
383  && !edgeIsCut_[findEdge(facei, v0, vMin1)]
384  )
385  {
386  cuts[cutI++] = vertToEVert(v0);
387  startFp = f.fcIndex(fp);
388  break;
389  }
390  }
391  }
392 
393  // Pass2: first edge cut not preceeded by point cut
394  if (startFp == -1)
395  {
396  forAll(f, fp)
397  {
398  label fp1 = f.fcIndex(fp);
399 
400  label v0 = f[fp];
401  label v1 = f[fp1];
402 
403  label edgeI = findEdge(facei, v0, v1);
404 
405  if (edgeIsCut_[edgeI] && !pointIsCut_[v0])
406  {
407  cuts[cutI++] = edgeToEVert(edgeI);
408  startFp = fp1;
409  break;
410  }
411  }
412  }
413 
414  // Pass3: nothing found so far. Either face is not cut at all or all
415  // vertices are cut. Start from 0.
416  if (startFp == -1)
417  {
418  startFp = 0;
419  }
420 
421  // Store all cuts starting from startFp;
422  label fp = startFp;
423 
424  bool allVerticesCut = true;
425 
426  forAll(f, i)
427  {
428  label fp1 = f.fcIndex(fp);
429 
430  // Get the three items: current vertex, next vertex and edge
431  // inbetween
432  label v0 = f[fp];
433  label v1 = f[fp1];
434  label edgeI = findEdge(facei, v0, v1);
435 
436  if (pointIsCut_[v0])
437  {
438  cuts[cutI++] = vertToEVert(v0);
439  }
440  else
441  {
442  // For check if all vertices have been cut (= illegal)
443  allVerticesCut = false;
444  }
445 
446  if (edgeIsCut_[edgeI])
447  {
448  cuts[cutI++] = edgeToEVert(edgeI);
449  }
450 
451  fp = fp1;
452  }
453 
454 
455  if (allVerticesCut)
456  {
458  << "Face " << facei << " vertices " << f
459  << " has all its vertices cut. Not cutting face." << endl;
460 
461  cutI = 0;
462  }
463 
464  // Remove duplicate starting point
465  if (cutI > 1 && cuts[cutI-1] == cuts[0])
466  {
467  cutI--;
468  }
469  cuts.setSize(cutI);
470  }
471 }
472 
473 
474 Foam::label Foam::cellCuts::findEdge
475 (
476  const label facei,
477  const label v0,
478  const label v1
479 ) const
480 {
481  const edgeList& edges = mesh().edges();
482 
483  const labelList& fEdges = mesh().faceEdges()[facei];
484 
485  forAll(fEdges, i)
486  {
487  const edge& e = edges[fEdges[i]];
488 
489  if
490  (
491  (e[0] == v0 && e[1] == v1)
492  || (e[0] == v1 && e[1] == v0)
493  )
494  {
495  return fEdges[i];
496  }
497  }
498 
499  return -1;
500 }
501 
502 
503 Foam::label Foam::cellCuts::loopFace
504 (
505  const label celli,
506  const labelList& loop
507 ) const
508 {
509  const cell& cFaces = mesh().cells()[celli];
510 
511  forAll(cFaces, cFacei)
512  {
513  label facei = cFaces[cFacei];
514 
515  const labelList& fEdges = mesh().faceEdges()[facei];
516  const face& f = mesh().faces()[facei];
517 
518  bool allOnFace = true;
519 
520  forAll(loop, i)
521  {
522  label cut = loop[i];
523 
524  if (isEdge(cut))
525  {
526  if (findIndex(fEdges, getEdge(cut)) == -1)
527  {
528  // Edge not on face. Skip face.
529  allOnFace = false;
530  break;
531  }
532  }
533  else
534  {
535  if (findIndex(f, getVertex(cut)) == -1)
536  {
537  // Vertex not on face. Skip face.
538  allOnFace = false;
539  break;
540  }
541  }
542  }
543 
544  if (allOnFace)
545  {
546  // Found face where all elements of loop are on the face.
547  return facei;
548  }
549  }
550  return -1;
551 }
552 
553 
554 bool Foam::cellCuts::walkPoint
555 (
556  const label celli,
557  const label startCut,
558 
559  const label exclude0,
560  const label exclude1,
561 
562  const label otherCut,
563 
564  label& nVisited,
565  labelList& visited
566 ) const
567 {
568  label vertI = getVertex(otherCut);
569 
570  const labelList& pFaces = mesh().pointFaces()[vertI];
571 
572  forAll(pFaces, pFacei)
573  {
574  label otherFacei = pFaces[pFacei];
575 
576  if
577  (
578  otherFacei != exclude0
579  && otherFacei != exclude1
580  && meshTools::faceOnCell(mesh(), celli, otherFacei)
581  )
582  {
583  label oldNVisited = nVisited;
584 
585  bool foundLoop =
586  walkCell
587  (
588  celli,
589  startCut,
590  otherFacei,
591  otherCut,
592  nVisited,
593  visited
594  );
595 
596  if (foundLoop)
597  {
598  return true;
599  }
600 
601  // No success. Restore state and continue
602  nVisited = oldNVisited;
603  }
604  }
605  return false;
606 }
607 
608 
609 bool Foam::cellCuts::crossEdge
610 (
611  const label celli,
612  const label startCut,
613  const label facei,
614  const label otherCut,
615 
616  label& nVisited,
617  labelList& visited
618 ) const
619 {
620  // Cross edge to other face
621  label edgeI = getEdge(otherCut);
622 
623  label otherFacei = meshTools::otherFace(mesh(), celli, facei, edgeI);
624 
625  // Store old state
626  label oldNVisited = nVisited;
627 
628  // Recurse to otherCut
629  bool foundLoop =
630  walkCell
631  (
632  celli,
633  startCut,
634  otherFacei,
635  otherCut,
636  nVisited,
637  visited
638  );
639 
640  if (foundLoop)
641  {
642  return true;
643  }
644  else
645  {
646  // No success. Restore state (i.e. backtrack)
647  nVisited = oldNVisited;
648 
649  return false;
650  }
651 }
652 
653 
654 bool Foam::cellCuts::addCut
655 (
656  const label celli,
657  const label cut,
658  label& nVisited,
659  labelList& visited
660 ) const
661 {
662  // Check for duplicate cuts.
663  if (findPartIndex(visited, nVisited, cut) != -1)
664  {
665  // Truncate (copy of) visited for ease of printing.
666  labelList truncVisited(visited);
667  truncVisited.setSize(nVisited);
668 
669  Pout<< "For cell " << celli << " : trying to add duplicate cut " << cut;
670  labelList cuts(1, cut);
671  writeCuts(Pout, cuts, loopWeights(cuts));
672 
673  Pout<< " to path:";
674  writeCuts(Pout, truncVisited, loopWeights(truncVisited));
675  Pout<< endl;
676 
677  return false;
678  }
679  else
680  {
681  visited[nVisited++] = cut;
682 
683  return true;
684  }
685 }
686 
687 
688 bool Foam::cellCuts::walkFace
689 (
690  const label celli,
691  const label startCut,
692  const label facei,
693  const label cut,
694 
695  label& lastCut,
696  label& beforeLastCut,
697  label& nVisited,
698  labelList& visited
699 ) const
700 {
701  const labelList& fCuts = faceCuts()[facei];
702 
703  if (fCuts.size() < 2)
704  {
705  return false;
706  }
707 
708  // Easy case : two cuts.
709  if (fCuts.size() == 2)
710  {
711  if (fCuts[0] == cut)
712  {
713  if (!addCut(celli, cut, nVisited, visited))
714  {
715  return false;
716  }
717 
718  beforeLastCut = cut;
719  lastCut = fCuts[1];
720 
721  return true;
722  }
723  else
724  {
725  if (!addCut(celli, cut, nVisited, visited))
726  {
727  return false;
728  }
729 
730  beforeLastCut = cut;
731  lastCut = fCuts[0];
732 
733  return true;
734  }
735  }
736 
737  // Harder case: more than 2 cuts on face.
738  // Should be start or end of string of cuts. Store all but last cut.
739  if (fCuts[0] == cut)
740  {
741  // Walk forward
742  for (label i = 0; i < fCuts.size()-1; i++)
743  {
744  if (!addCut(celli, fCuts[i], nVisited, visited))
745  {
746  return false;
747  }
748  }
749  beforeLastCut = fCuts[fCuts.size()-2];
750  lastCut = fCuts[fCuts.size()-1];
751  }
752  else if (fCuts[fCuts.size()-1] == cut)
753  {
754  for (label i = fCuts.size()-1; i >= 1; --i)
755  {
756  if (!addCut(celli, fCuts[i], nVisited, visited))
757  {
758  return false;
759  }
760  }
761  beforeLastCut = fCuts[1];
762  lastCut = fCuts[0];
763  }
764  else
765  {
767  << "In middle of cut. cell:" << celli << " face:" << facei
768  << " cuts:" << fCuts << " current cut:" << cut << endl;
769 
770  return false;
771  }
772 
773  return true;
774 }
775 
776 
777 
778 bool Foam::cellCuts::walkCell
779 (
780  const label celli,
781  const label startCut,
782  const label facei,
783  const label cut,
784 
785  label& nVisited,
786  labelList& visited
787 ) const
788 {
789  // Walk across face, storing cuts. Return last two cuts visited.
790  label lastCut = -1;
791  label beforeLastCut = -1;
792 
793 
794  if (debug & 2)
795  {
796  Pout<< "For cell:" << celli << " walked across face " << facei
797  << " from cut ";
798  labelList cuts(1, cut);
799  writeCuts(Pout, cuts, loopWeights(cuts));
800  Pout<< endl;
801  }
802 
803  bool validWalk = walkFace
804  (
805  celli,
806  startCut,
807  facei,
808  cut,
809 
810  lastCut,
811  beforeLastCut,
812  nVisited,
813  visited
814  );
815 
816  if (!validWalk)
817  {
818  return false;
819  }
820 
821  if (debug & 2)
822  {
823  Pout<< " to last cut ";
824  labelList cuts(1, lastCut);
825  writeCuts(Pout, cuts, loopWeights(cuts));
826  Pout<< endl;
827  }
828 
829  // Check if starting point reached.
830  if (lastCut == startCut)
831  {
832  if (nVisited >= 3)
833  {
834  if (debug & 2)
835  {
836  // Truncate (copy of) visited for ease of printing.
837  labelList truncVisited(visited);
838  truncVisited.setSize(nVisited);
839 
840  Pout<< "For cell " << celli << " : found closed path:";
841  writeCuts(Pout, truncVisited, loopWeights(truncVisited));
842  Pout<< " closed by " << lastCut << endl;
843  }
844 
845  return true;
846  }
847  else
848  {
849  // Loop but too short. Probably folded back on itself.
850  return false;
851  }
852  }
853 
854 
855  // Check last cut and one before that to determine type.
856 
857  // From beforeLastCut to lastCut:
858  // - from edge to edge
859  // (- always not along existing edge)
860  // - from edge to vertex
861  // - not along existing edge
862  // - along existing edge: not allowed?
863  // - from vertex to edge
864  // - not along existing edge
865  // - along existing edge. Not allowed. See above.
866  // - from vertex to vertex
867  // - not along existing edge
868  // - along existing edge
869 
870  if (isEdge(beforeLastCut))
871  {
872  if (isEdge(lastCut))
873  {
874  // beforeLastCut=edge, lastCut=edge.
875 
876  // Cross lastCut (=edge) into face not facei
877  return crossEdge
878  (
879  celli,
880  startCut,
881  facei,
882  lastCut,
883  nVisited,
884  visited
885  );
886  }
887  else
888  {
889  // beforeLastCut=edge, lastCut=vertex.
890 
891  // Go from lastCut to all connected faces (except facei)
892  return walkPoint
893  (
894  celli,
895  startCut,
896  facei, // exclude0: don't cross back on current face
897  -1, // exclude1
898  lastCut,
899  nVisited,
900  visited
901  );
902  }
903  }
904  else
905  {
906  if (isEdge(lastCut))
907  {
908  // beforeLastCut=vertex, lastCut=edge.
909  return crossEdge
910  (
911  celli,
912  startCut,
913  facei,
914  lastCut,
915  nVisited,
916  visited
917  );
918  }
919  else
920  {
921  // beforeLastCut=vertex, lastCut=vertex. Check if along existing
922  // edge.
923  label edgeI =
924  findEdge
925  (
926  facei,
927  getVertex(beforeLastCut),
928  getVertex(lastCut)
929  );
930 
931  if (edgeI != -1)
932  {
933  // Cut along existing edge. So is in fact on two faces.
934  // Get faces on both sides of the edge to make
935  // sure we dont fold back on to those.
936 
937  label f0, f1;
938  meshTools::getEdgeFaces(mesh(), celli, edgeI, f0, f1);
939 
940  return walkPoint
941  (
942  celli,
943  startCut,
944  f0,
945  f1,
946  lastCut,
947  nVisited,
948  visited
949  );
950 
951  }
952  else
953  {
954  // Cross cut across face.
955  return walkPoint
956  (
957  celli,
958  startCut,
959  facei, // face to exclude
960  -1, // face to exclude
961  lastCut,
962  nVisited,
963  visited
964  );
965  }
966  }
967  }
968 }
969 
970 
971 void Foam::cellCuts::calcCellLoops(const labelList& cutCells)
972 {
973  // Determine for every cut cell the loop (= face) it is cut by. Done by
974  // starting from a cut edge or cut vertex and walking across faces, from
975  // cut to cut, until starting cut hit.
976  // If multiple loops are possible across a cell circumference takes the
977  // first one found.
978 
979  // Calculate cuts per face.
980  const labelListList& allFaceCuts = faceCuts();
981 
982  // Per cell the number of faces with valid cuts. Is used as quick
983  // rejection to see if cell can be cut.
984  labelList nCutFaces(mesh().nCells(), 0);
985 
986  forAll(allFaceCuts, facei)
987  {
988  const labelList& fCuts = allFaceCuts[facei];
989 
990  if (fCuts.size() == mesh().faces()[facei].size())
991  {
992  // Too many cuts on face. WalkCell would get very upset so disable.
993  nCutFaces[mesh().faceOwner()[facei]] = labelMin;
994 
995  if (mesh().isInternalFace(facei))
996  {
997  nCutFaces[mesh().faceNeighbour()[facei]] = labelMin;
998  }
999  }
1000  else if (fCuts.size() >= 2)
1001  {
1002  // Could be valid cut. Update count for owner and neighbour.
1003  nCutFaces[mesh().faceOwner()[facei]]++;
1004 
1005  if (mesh().isInternalFace(facei))
1006  {
1007  nCutFaces[mesh().faceNeighbour()[facei]]++;
1008  }
1009  }
1010  }
1011 
1012 
1013  // Stack of visited cuts (nVisited used as stack pointer)
1014  // Size big enough.
1015  labelList visited(mesh().nPoints());
1016 
1017  forAll(cutCells, i)
1018  {
1019  label celli = cutCells[i];
1020 
1021  bool validLoop = false;
1022 
1023  // Quick rejection: has enough faces that are cut?
1024  if (nCutFaces[celli] >= 3)
1025  {
1026  const labelList& cFaces = mesh().cells()[celli];
1027 
1028  if (debug & 2)
1029  {
1030  Pout<< "cell:" << celli << " cut faces:" << endl;
1031  forAll(cFaces, i)
1032  {
1033  label facei = cFaces[i];
1034  const labelList& fCuts = allFaceCuts[facei];
1035 
1036  Pout<< " face:" << facei << " cuts:";
1037  writeCuts(Pout, fCuts, loopWeights(fCuts));
1038  Pout<< endl;
1039  }
1040  }
1041 
1042  label nVisited = 0;
1043 
1044  // Determine the first cut face to start walking from.
1045  forAll(cFaces, cFacei)
1046  {
1047  label facei = cFaces[cFacei];
1048 
1049  const labelList& fCuts = allFaceCuts[facei];
1050 
1051  // Take first or last cut of multiple on face.
1052  // Note that in calcFaceCuts
1053  // we have already made sure this is the start or end of a
1054  // string of cuts.
1055  if (fCuts.size() >= 2)
1056  {
1057  // Try walking from start of fCuts.
1058  nVisited = 0;
1059 
1060  if (debug & 2)
1061  {
1062  Pout<< "cell:" << celli
1063  << " start walk at face:" << facei
1064  << " cut:";
1065  labelList cuts(1, fCuts[0]);
1066  writeCuts(Pout, cuts, loopWeights(cuts));
1067  Pout<< endl;
1068  }
1069 
1070  validLoop =
1071  walkCell
1072  (
1073  celli,
1074  fCuts[0],
1075  facei,
1076  fCuts[0],
1077 
1078  nVisited,
1079  visited
1080  );
1081 
1082  if (validLoop)
1083  {
1084  break;
1085  }
1086 
1087  // No need to try and walk from end of cuts since
1088  // all paths already tried by walkCell.
1089  }
1090  }
1091 
1092  if (validLoop)
1093  {
1094  // Copy nVisited elements out of visited (since visited is
1095  // never truncated for efficiency reasons)
1096 
1097  labelList& loop = cellLoops_[celli];
1098 
1099  loop.setSize(nVisited);
1100 
1101  for (label i = 0; i < nVisited; i++)
1102  {
1103  loop[i] = visited[i];
1104  }
1105  }
1106  else
1107  {
1108  // Invalid loop. Leave cellLoops_[celli] zero size which
1109  // flags this.
1110  Pout<< "calcCellLoops(const labelList&) : did not find valid"
1111  << " loop for cell " << celli << endl;
1112  // Dump cell and cuts on cell.
1113  writeUncutOBJ(".", celli);
1114 
1115  cellLoops_[celli].setSize(0);
1116  }
1117  }
1118  else
1119  {
1120  //Pout<< "calcCellLoops(const labelList&) : did not find valid"
1121  // << " loop for cell " << celli << " since not enough cut faces"
1122  // << endl;
1123  cellLoops_[celli].setSize(0);
1124  }
1125  }
1126 }
1127 
1128 
1129 void Foam::cellCuts::walkEdges
1130 (
1131  const label celli,
1132  const label pointi,
1133  const label status,
1134 
1135  Map<label>& edgeStatus,
1136  Map<label>& pointStatus
1137 ) const
1138 {
1139  if (pointStatus.insert(pointi, status))
1140  {
1141  // First visit to pointi
1142 
1143  const labelList& pEdges = mesh().pointEdges()[pointi];
1144 
1145  forAll(pEdges, pEdgeI)
1146  {
1147  label edgeI = pEdges[pEdgeI];
1148 
1149  if
1150  (
1151  meshTools::edgeOnCell(mesh(), celli, edgeI)
1152  && edgeStatus.insert(edgeI, status)
1153  )
1154  {
1155  // First visit to edgeI so recurse.
1156 
1157  label v2 = mesh().edges()[edgeI].otherVertex(pointi);
1158 
1159  walkEdges(celli, v2, status, edgeStatus, pointStatus);
1160  }
1161  }
1162  }
1163 }
1164 
1165 
1168  const labelList& cellPoints,
1169  const labelList& anchorPoints,
1170  const labelList& loop
1171 ) const
1172 {
1173  labelList newElems(cellPoints.size());
1174  label newElemI = 0;
1175 
1176  forAll(cellPoints, i)
1177  {
1178  label pointi = cellPoints[i];
1179 
1180  if
1181  (
1182  findIndex(anchorPoints, pointi) == -1
1183  && findIndex(loop, vertToEVert(pointi)) == -1
1184  )
1185  {
1186  newElems[newElemI++] = pointi;
1187  }
1188  }
1189 
1190  newElems.setSize(newElemI);
1191 
1192  return newElems;
1193 }
1194 
1195 
1196 bool Foam::cellCuts::loopAnchorConsistent
1197 (
1198  const label celli,
1199  const pointField& loopPts,
1200  const labelList& anchorPoints
1201 ) const
1202 {
1203  // Create identity face for ease of calculation of normal etc.
1204  face f(identity(loopPts.size()));
1205 
1206  vector normal = f.normal(loopPts);
1207  point ctr = f.centre(loopPts);
1208 
1209 
1210  // Get average position of anchor points.
1211  vector avg(Zero);
1212 
1213  forAll(anchorPoints, ptI)
1214  {
1215  avg += mesh().points()[anchorPoints[ptI]];
1216  }
1217  avg /= anchorPoints.size();
1218 
1219 
1220  if (((avg - ctr) & normal) > 0)
1221  {
1222  return true;
1223  }
1224  else
1225  {
1226  return false;
1227  }
1228 }
1229 
1230 
1231 bool Foam::cellCuts::calcAnchors
1232 (
1233  const label celli,
1234  const labelList& loop,
1235  const pointField& loopPts,
1236 
1237  labelList& anchorPoints
1238 ) const
1239 {
1240  const edgeList& edges = mesh().edges();
1241 
1242  const labelList& cPoints = mesh().cellPoints()[celli];
1243  const labelList& cEdges = mesh().cellEdges()[celli];
1244  const cell& cFaces = mesh().cells()[celli];
1245 
1246  // Points on loop
1247 
1248  // Status of point:
1249  // 0 - on loop
1250  // 1 - point set 1
1251  // 2 - point set 2
1252  Map<label> pointStatus(2*cPoints.size());
1253  Map<label> edgeStatus(2*cEdges.size());
1254 
1255  // Mark loop vertices
1256  forAll(loop, i)
1257  {
1258  label cut = loop[i];
1259 
1260  if (isEdge(cut))
1261  {
1262  edgeStatus.insert(getEdge(cut), 0);
1263  }
1264  else
1265  {
1266  pointStatus.insert(getVertex(cut), 0);
1267  }
1268  }
1269  // Since edges between two cut vertices have not been marked, mark them
1270  // explicitly
1271  forAll(cEdges, i)
1272  {
1273  label edgeI = cEdges[i];
1274  const edge& e = edges[edgeI];
1275 
1276  if (pointStatus.found(e[0]) && pointStatus.found(e[1]))
1277  {
1278  edgeStatus.insert(edgeI, 0);
1279  }
1280  }
1281 
1282 
1283  // Find uncut starting vertex
1284  label uncutIndex = firstUnique(cPoints, pointStatus);
1285 
1286  if (uncutIndex == -1)
1287  {
1289  << "Invalid loop " << loop << " for cell " << celli << endl
1290  << "Can not find point on cell which is not cut by loop."
1291  << endl;
1292 
1293  writeOBJ(".", celli, loopPts, labelList(0));
1294 
1295  return false;
1296  }
1297 
1298  // Walk unset vertices and edges and mark with 1 in pointStatus, edgeStatus
1299  walkEdges(celli, cPoints[uncutIndex], 1, edgeStatus, pointStatus);
1300 
1301  // Find new uncut starting vertex
1302  uncutIndex = firstUnique(cPoints, pointStatus);
1303 
1304  if (uncutIndex == -1)
1305  {
1306  // All vertices either in loop or in anchor. So split is along single
1307  // face.
1309  << "Invalid loop " << loop << " for cell " << celli << endl
1310  << "All vertices of cell are either in loop or in anchor set"
1311  << endl;
1312 
1313  writeOBJ(".", celli, loopPts, labelList(0));
1314 
1315  return false;
1316  }
1317 
1318  // Walk unset vertices and edges and mark with 2. These are the
1319  // pointset 2.
1320  walkEdges(celli, cPoints[uncutIndex], 2, edgeStatus, pointStatus);
1321 
1322  // Collect both sets in lists.
1323  DynamicList<label> connectedPoints(cPoints.size());
1324  DynamicList<label> otherPoints(cPoints.size());
1325 
1326  forAllConstIter(Map<label>, pointStatus, iter)
1327  {
1328  if (iter() == 1)
1329  {
1330  connectedPoints.append(iter.key());
1331  }
1332  else if (iter() == 2)
1333  {
1334  otherPoints.append(iter.key());
1335  }
1336  }
1337  connectedPoints.shrink();
1338  otherPoints.shrink();
1339 
1340  // Check that all points have been used.
1341  uncutIndex = firstUnique(cPoints, pointStatus);
1342 
1343  if (uncutIndex != -1)
1344  {
1346  << "Invalid loop " << loop << " for cell " << celli
1347  << " since it splits the cell into more than two cells" << endl;
1348 
1349  writeOBJ(".", celli, loopPts, connectedPoints);
1350 
1351  return false;
1352  }
1353 
1354 
1355  // Check that both parts (connectedPoints, otherPoints) have enough faces.
1356  labelHashSet connectedFaces(2*cFaces.size());
1357  labelHashSet otherFaces(2*cFaces.size());
1358 
1359  forAllConstIter(Map<label>, pointStatus, iter)
1360  {
1361  label pointi = iter.key();
1362 
1363  const labelList& pFaces = mesh().pointFaces()[pointi];
1364 
1365  if (iter() == 1)
1366  {
1367  forAll(pFaces, pFacei)
1368  {
1369  if (meshTools::faceOnCell(mesh(), celli, pFaces[pFacei]))
1370  {
1371  connectedFaces.insert(pFaces[pFacei]);
1372  }
1373  }
1374  }
1375  else if (iter() == 2)
1376  {
1377  forAll(pFaces, pFacei)
1378  {
1379  if (meshTools::faceOnCell(mesh(), celli, pFaces[pFacei]))
1380  {
1381  otherFaces.insert(pFaces[pFacei]);
1382  }
1383  }
1384  }
1385  }
1386 
1387  if (connectedFaces.size() < 3)
1388  {
1390  << "Invalid loop " << loop << " for cell " << celli
1391  << " since would have too few faces on one side." << nl
1392  << "All faces:" << cFaces << endl;
1393 
1394  writeOBJ(".", celli, loopPts, connectedPoints);
1395 
1396  return false;
1397  }
1398 
1399  if (otherFaces.size() < 3)
1400  {
1402  << "Invalid loop " << loop << " for cell " << celli
1403  << " since would have too few faces on one side." << nl
1404  << "All faces:" << cFaces << endl;
1405 
1406  writeOBJ(".", celli, loopPts, otherPoints);
1407 
1408  return false;
1409  }
1410 
1411 
1412  // Check that faces are split into two regions and not more.
1413  // When walking across the face the only transition of pointStatus is
1414  // from set1 to loop to set2 or back. Not allowed is from set1 to loop to
1415  // set1.
1416  {
1417  forAll(cFaces, i)
1418  {
1419  label facei = cFaces[i];
1420 
1421  const face& f = mesh().faces()[facei];
1422 
1423  bool hasSet1 = false;
1424  bool hasSet2 = false;
1425 
1426  label prevStat = pointStatus[f[0]];
1427 
1428  forAll(f, fp)
1429  {
1430  label v0 = f[fp];
1431  label pStat = pointStatus[v0];
1432 
1433  if (pStat == prevStat)
1434  {
1435  }
1436  else if (pStat == 0)
1437  {
1438  // Loop.
1439  }
1440  else if (pStat == 1)
1441  {
1442  if (hasSet1)
1443  {
1444  // Second occurence of set1.
1446  << "Invalid loop " << loop << " for cell " << celli
1447  << " since face " << f << " would be split into"
1448  << " more than two faces" << endl;
1449 
1450  writeOBJ(".", celli, loopPts, otherPoints);
1451 
1452  return false;
1453  }
1454 
1455  hasSet1 = true;
1456  }
1457  else if (pStat == 2)
1458  {
1459  if (hasSet2)
1460  {
1461  // Second occurence of set1.
1463  << "Invalid loop " << loop << " for cell " << celli
1464  << " since face " << f << " would be split into"
1465  << " more than two faces" << endl;
1466 
1467  writeOBJ(".", celli, loopPts, otherPoints);
1468 
1469  return false;
1470  }
1471 
1472  hasSet2 = true;
1473  }
1474  else
1475  {
1477  << abort(FatalError);
1478  }
1479 
1480  prevStat = pStat;
1481 
1482 
1483  label v1 = f.nextLabel(fp);
1484  label edgeI = findEdge(facei, v0, v1);
1485 
1486  label eStat = edgeStatus[edgeI];
1487 
1488  if (eStat == prevStat)
1489  {
1490  }
1491  else if (eStat == 0)
1492  {
1493  // Loop.
1494  }
1495  else if (eStat == 1)
1496  {
1497  if (hasSet1)
1498  {
1499  // Second occurence of set1.
1501  << "Invalid loop " << loop << " for cell " << celli
1502  << " since face " << f << " would be split into"
1503  << " more than two faces" << endl;
1504 
1505  writeOBJ(".", celli, loopPts, otherPoints);
1506 
1507  return false;
1508  }
1509 
1510  hasSet1 = true;
1511  }
1512  else if (eStat == 2)
1513  {
1514  if (hasSet2)
1515  {
1516  // Second occurence of set1.
1518  << "Invalid loop " << loop << " for cell " << celli
1519  << " since face " << f << " would be split into"
1520  << " more than two faces" << endl;
1521 
1522  writeOBJ(".", celli, loopPts, otherPoints);
1523 
1524  return false;
1525  }
1526 
1527  hasSet2 = true;
1528  }
1529  prevStat = eStat;
1530  }
1531  }
1532  }
1533 
1534 
1535 
1536 
1537  // Check which one of point sets to use.
1538  bool loopOk = loopAnchorConsistent(celli, loopPts, connectedPoints);
1539 
1540  //if (debug)
1541  {
1542  // Additional check: are non-anchor points on other side?
1543  bool otherLoopOk = loopAnchorConsistent(celli, loopPts, otherPoints);
1544 
1545  if (loopOk == otherLoopOk)
1546  {
1547  // Both sets of points are supposedly on the same side as the
1548  // loop normal. Oops.
1549 
1551  << " For cell:" << celli
1552  << " achorpoints and nonanchorpoints are geometrically"
1553  << " on same side!" << endl
1554  << "cellPoints:" << cPoints << endl
1555  << "loop:" << loop << endl
1556  << "anchors:" << connectedPoints << endl
1557  << "otherPoints:" << otherPoints << endl;
1558 
1559  writeOBJ(".", celli, loopPts, connectedPoints);
1560  }
1561  }
1562 
1563  if (loopOk)
1564  {
1565  // connectedPoints on 'outside' of loop so these are anchor points
1566  anchorPoints.transfer(connectedPoints);
1567  connectedPoints.clear();
1568  }
1569  else
1570  {
1571  anchorPoints.transfer(otherPoints);
1572  otherPoints.clear();
1573  }
1574  return true;
1575 }
1576 
1577 
1578 Foam::pointField Foam::cellCuts::loopPoints
1579 (
1580  const labelList& loop,
1581  const scalarField& loopWeights
1582 ) const
1583 {
1584  pointField loopPts(loop.size());
1585 
1586  forAll(loop, fp)
1587  {
1588  loopPts[fp] = coord(loop[fp], loopWeights[fp]);
1589  }
1590  return loopPts;
1591 }
1592 
1593 
1594 Foam::scalarField Foam::cellCuts::loopWeights(const labelList& loop) const
1595 {
1596  scalarField weights(loop.size());
1597 
1598  forAll(loop, fp)
1599  {
1600  label cut = loop[fp];
1601 
1602  if (isEdge(cut))
1603  {
1604  label edgeI = getEdge(cut);
1605 
1606  weights[fp] = edgeWeight_[edgeI];
1607  }
1608  else
1609  {
1610  weights[fp] = -GREAT;
1611  }
1612  }
1613  return weights;
1614 }
1615 
1616 
1617 bool Foam::cellCuts::validEdgeLoop
1618 (
1619  const labelList& loop,
1620  const scalarField& loopWeights
1621 ) const
1622 {
1623  forAll(loop, fp)
1624  {
1625  label cut = loop[fp];
1626 
1627  if (isEdge(cut))
1628  {
1629  label edgeI = getEdge(cut);
1630 
1631  // Check: cut compatible only if can be snapped to existing one.
1632  if (edgeIsCut_[edgeI])
1633  {
1634  scalar edgeLen = mesh().edges()[edgeI].mag(mesh().points());
1635 
1636  if
1637  (
1638  mag(loopWeights[fp] - edgeWeight_[edgeI])
1639  > geomCellLooper::snapTol()*edgeLen
1640  )
1641  {
1642  // Positions differ too much->would create two cuts.
1643  return false;
1644  }
1645  }
1646  }
1647  }
1648  return true;
1649 }
1650 
1651 
1652 Foam::label Foam::cellCuts::countFaceCuts
1653 (
1654  const label facei,
1655  const labelList& loop
1656 ) const
1657 {
1658  // Includes cuts through vertices and through edges.
1659  // Assumes that if edge is cut both in edgeIsCut and in loop that the
1660  // position of the cut is the same.
1661 
1662  label nCuts = 0;
1663 
1664  // Count cut vertices
1665  const face& f = mesh().faces()[facei];
1666 
1667  forAll(f, fp)
1668  {
1669  label vertI = f[fp];
1670 
1671  // Vertex already cut or mentioned in current loop.
1672  if
1673  (
1674  pointIsCut_[vertI]
1675  || (findIndex(loop, vertToEVert(vertI)) != -1)
1676  )
1677  {
1678  nCuts++;
1679  }
1680  }
1681 
1682  // Count cut edges.
1683  const labelList& fEdges = mesh().faceEdges()[facei];
1684 
1685  forAll(fEdges, fEdgeI)
1686  {
1687  label edgeI = fEdges[fEdgeI];
1688 
1689  // Edge already cut or mentioned in current loop.
1690  if
1691  (
1692  edgeIsCut_[edgeI]
1693  || (findIndex(loop, edgeToEVert(edgeI)) != -1)
1694  )
1695  {
1696  nCuts++;
1697  }
1698  }
1699 
1700  return nCuts;
1701 }
1702 
1703 
1704 bool Foam::cellCuts::conservativeValidLoop
1705 (
1706  const label celli,
1707  const labelList& loop
1708 ) const
1709 {
1710 
1711  if (loop.size() < 2)
1712  {
1713  return false;
1714  }
1715 
1716  forAll(loop, cutI)
1717  {
1718  if (isEdge(loop[cutI]))
1719  {
1720  label edgeI = getEdge(loop[cutI]);
1721 
1722  if (edgeIsCut_[edgeI])
1723  {
1724  // edge compatibility already checked.
1725  }
1726  else
1727  {
1728  // Quick rejection: vertices of edge itself cannot be cut.
1729  const edge& e = mesh().edges()[edgeI];
1730 
1731  if (pointIsCut_[e.start()] || pointIsCut_[e.end()])
1732  {
1733  return false;
1734  }
1735 
1736 
1737  // Check faces using this edge
1738  const labelList& eFaces = mesh().edgeFaces()[edgeI];
1739 
1740  forAll(eFaces, eFacei)
1741  {
1742  label nCuts = countFaceCuts(eFaces[eFacei], loop);
1743 
1744  if (nCuts > 2)
1745  {
1746  return false;
1747  }
1748  }
1749  }
1750  }
1751  else
1752  {
1753  // Vertex cut
1754 
1755  label vertI = getVertex(loop[cutI]);
1756 
1757  if (!pointIsCut_[vertI])
1758  {
1759  // New cut through vertex.
1760 
1761  // Check edges using vertex.
1762  const labelList& pEdges = mesh().pointEdges()[vertI];
1763 
1764  forAll(pEdges, pEdgeI)
1765  {
1766  label edgeI = pEdges[pEdgeI];
1767 
1768  if (edgeIsCut_[edgeI])
1769  {
1770  return false;
1771  }
1772  }
1773 
1774  // Check faces using vertex.
1775  const labelList& pFaces = mesh().pointFaces()[vertI];
1776 
1777  forAll(pFaces, pFacei)
1778  {
1779  label nCuts = countFaceCuts(pFaces[pFacei], loop);
1780 
1781  if (nCuts > 2)
1782  {
1783  return false;
1784  }
1785  }
1786  }
1787  }
1788  }
1789 
1790  // All ok.
1791  return true;
1792 }
1793 
1794 
1795 bool Foam::cellCuts::validLoop
1796 (
1797  const label celli,
1798  const labelList& loop,
1799  const scalarField& loopWeights,
1800 
1801  Map<edge>& newFaceSplitCut,
1802  labelList& anchorPoints
1803 ) const
1804 {
1805  // Determine compatibility of loop with existing cut pattern. Does not use
1806  // derived cut-addressing (faceCuts), only pointIsCut, edgeIsCut.
1807  // Adds any cross-cuts found to newFaceSplitCut and sets cell points on
1808  // one side of the loop in anchorPoints.
1809 
1810  if (loop.size() < 2)
1811  {
1812  return false;
1813  }
1814 
1815  if (debug & 4)
1816  {
1817  // Allow as fallback the 'old' loop checking where only a single
1818  // cut per face is allowed.
1819  if (!conservativeValidLoop(celli, loop))
1820  {
1821  Info << "Invalid conservative loop: " << loop << endl;
1822  return false;
1823  }
1824  }
1825 
1826  forAll(loop, fp)
1827  {
1828  label cut = loop[fp];
1829  label nextCut = loop[(fp+1) % loop.size()];
1830 
1831  // Label (if any) of face cut (so cut not along existing edge)
1832  label meshFacei = -1;
1833 
1834  if (isEdge(cut))
1835  {
1836  label edgeI = getEdge(cut);
1837 
1838  // Look one cut ahead to find if it is along existing edge.
1839 
1840  if (isEdge(nextCut))
1841  {
1842  // From edge to edge -> cross cut
1843  label nextEdgeI = getEdge(nextCut);
1844 
1845  // Find face and mark as to be split.
1846  meshFacei = edgeEdgeToFace(celli, edgeI, nextEdgeI);
1847 
1848  if (meshFacei == -1)
1849  {
1850  // Can't find face using both cut edges.
1851  return false;
1852  }
1853  }
1854  else
1855  {
1856  // From edge to vertex -> cross cut only if no existing edge.
1857 
1858  label nextVertI = getVertex(nextCut);
1859  const edge& e = mesh().edges()[edgeI];
1860 
1861  if (e.start() != nextVertI && e.end() != nextVertI)
1862  {
1863  // New edge. Find face and mark as to be split.
1864  meshFacei = edgeVertexToFace(celli, edgeI, nextVertI);
1865 
1866  if (meshFacei == -1)
1867  {
1868  // Can't find face. Ilegal.
1869  return false;
1870  }
1871  }
1872  }
1873  }
1874  else
1875  {
1876  // Cut is vertex.
1877  label vertI = getVertex(cut);
1878 
1879  if (isEdge(nextCut))
1880  {
1881  // From vertex to edge -> cross cut only if no existing edge.
1882  label nextEdgeI = getEdge(nextCut);
1883 
1884  const edge& nextE = mesh().edges()[nextEdgeI];
1885 
1886  if (nextE.start() != vertI && nextE.end() != vertI)
1887  {
1888  // Should be cross cut. Find face.
1889  meshFacei = edgeVertexToFace(celli, nextEdgeI, vertI);
1890 
1891  if (meshFacei == -1)
1892  {
1893  return false;
1894  }
1895  }
1896  }
1897  else
1898  {
1899  // From vertex to vertex -> cross cut only if no existing edge.
1900  label nextVertI = getVertex(nextCut);
1901 
1902  if (meshTools::findEdge(mesh(), vertI, nextVertI) == -1)
1903  {
1904  // New edge. Find face.
1905  meshFacei = vertexVertexToFace(celli, vertI, nextVertI);
1906 
1907  if (meshFacei == -1)
1908  {
1909  return false;
1910  }
1911  }
1912  }
1913  }
1914 
1915  if (meshFacei != -1)
1916  {
1917  // meshFacei is cut across along cut-nextCut (so not along existing
1918  // edge). Check if this is compatible with existing pattern.
1919  edge cutEdge(cut, nextCut);
1920 
1921  Map<edge>::const_iterator iter = faceSplitCut_.find(meshFacei);
1922 
1923  if (iter == faceSplitCut_.end())
1924  {
1925  // Face not yet cut so insert.
1926  newFaceSplitCut.insert(meshFacei, cutEdge);
1927  }
1928  else
1929  {
1930  // Face already cut. Ok if same edge.
1931  if (iter() != cutEdge)
1932  {
1933  return false;
1934  }
1935  }
1936  }
1937  }
1938 
1939  // Is there a face on which all cuts are?
1940  label faceContainingLoop = loopFace(celli, loop);
1941 
1942  if (faceContainingLoop != -1)
1943  {
1945  << "Found loop on cell " << celli << " with all points"
1946  << " on face " << faceContainingLoop << endl;
1947 
1948  //writeOBJ(".", celli, loopPoints(loop, loopWeights), labelList(0));
1949 
1950  return false;
1951  }
1952 
1953  // Calculate anchor points
1954  // Final success is determined by whether anchor points can be determined.
1955  return calcAnchors
1956  (
1957  celli,
1958  loop,
1959  loopPoints(loop, loopWeights),
1960  anchorPoints
1961  );
1962 }
1963 
1964 
1965 void Foam::cellCuts::setFromCellLoops()
1966 {
1967  // 'Uncut' edges/vertices that are not used in loops.
1968  pointIsCut_ = false;
1969 
1970  edgeIsCut_ = false;
1971 
1972  faceSplitCut_.clear();
1973 
1974  forAll(cellLoops_, celli)
1975  {
1976  const labelList& loop = cellLoops_[celli];
1977 
1978  if (loop.size())
1979  {
1980  // Storage for cross-face cuts
1981  Map<edge> faceSplitCuts(loop.size());
1982  // Storage for points on one side of cell.
1983  labelList anchorPoints;
1984 
1985  if
1986  (
1987  !validLoop
1988  (
1989  celli,
1990  loop,
1991  loopWeights(loop),
1992  faceSplitCuts,
1993  anchorPoints
1994  )
1995  )
1996  {
1997  //writeOBJ(".", celli, loopPoints(celli), anchorPoints);
1998 
2000  << "Illegal loop " << loop
2001  << " when recreating cut-addressing"
2002  << " from existing cellLoops for cell " << celli
2003  << endl;
2004 
2005  cellLoops_[celli].setSize(0);
2006  cellAnchorPoints_[celli].setSize(0);
2007  }
2008  else
2009  {
2010  // Copy anchor points.
2011  cellAnchorPoints_[celli].transfer(anchorPoints);
2012 
2013  // Copy faceSplitCuts into overall faceSplit info.
2014  forAllConstIter(Map<edge>, faceSplitCuts, iter)
2015  {
2016  faceSplitCut_.insert(iter.key(), iter());
2017  }
2018 
2019  // Update edgeIsCut, pointIsCut information
2020  forAll(loop, cutI)
2021  {
2022  label cut = loop[cutI];
2023 
2024  if (isEdge(cut))
2025  {
2026  edgeIsCut_[getEdge(cut)] = true;
2027  }
2028  else
2029  {
2030  pointIsCut_[getVertex(cut)] = true;
2031  }
2032  }
2033  }
2034  }
2035  }
2036 
2037  // Reset edge weights
2038  forAll(edgeIsCut_, edgeI)
2039  {
2040  if (!edgeIsCut_[edgeI])
2041  {
2042  // Weight not used. Set to illegal value to make any use fall over.
2043  edgeWeight_[edgeI] = -GREAT;
2044  }
2045  }
2046 }
2047 
2048 
2049 bool Foam::cellCuts::setFromCellLoop
2050 (
2051  const label celli,
2052  const labelList& loop,
2053  const scalarField& loopWeights
2054 )
2055 {
2056  // Update basic cut information from single cellLoop. Returns true if loop
2057  // was valid.
2058 
2059  // Dump loop for debugging.
2060  if (debug)
2061  {
2062  OFstream str("last_cell.obj");
2063 
2064  str<< "# edges of cell " << celli << nl;
2065 
2067  (
2068  str,
2069  mesh().cells(),
2070  mesh().faces(),
2071  mesh().points(),
2072  labelList(1, celli)
2073  );
2074 
2075 
2076  OFstream loopStr("last_loop.obj");
2077 
2078  loopStr<< "# looppoints for cell " << celli << nl;
2079 
2080  pointField pointsOfLoop = loopPoints(loop, loopWeights);
2081 
2082  forAll(pointsOfLoop, i)
2083  {
2084  meshTools::writeOBJ(loopStr, pointsOfLoop[i]);
2085  }
2086 
2087  str << 'l';
2088 
2089  forAll(pointsOfLoop, i)
2090  {
2091  str << ' ' << i + 1;
2092  }
2093  str << ' ' << 1 << nl;
2094  }
2095 
2096  bool okLoop = false;
2097 
2098  if (validEdgeLoop(loop, loopWeights))
2099  {
2100  // Storage for cuts across face
2101  Map<edge> faceSplitCuts(loop.size());
2102  // Storage for points on one side of cell
2103  labelList anchorPoints;
2104 
2105  okLoop =
2106  validLoop(celli, loop, loopWeights, faceSplitCuts, anchorPoints);
2107 
2108  if (okLoop)
2109  {
2110  // Valid loop. Copy cellLoops and anchorPoints
2111  cellLoops_[celli] = loop;
2112  cellAnchorPoints_[celli].transfer(anchorPoints);
2113 
2114  // Copy split cuts
2115  forAllConstIter(Map<edge>, faceSplitCuts, iter)
2116  {
2117  faceSplitCut_.insert(iter.key(), iter());
2118  }
2119 
2120 
2121  // Update edgeIsCut, pointIsCut information
2122  forAll(loop, cutI)
2123  {
2124  label cut = loop[cutI];
2125 
2126  if (isEdge(cut))
2127  {
2128  label edgeI = getEdge(cut);
2129 
2130  edgeIsCut_[edgeI] = true;
2131 
2132  edgeWeight_[edgeI] = loopWeights[cutI];
2133  }
2134  else
2135  {
2136  label vertI = getVertex(cut);
2137 
2138  pointIsCut_[vertI] = true;
2139  }
2140  }
2141  }
2142  }
2143 
2144  return okLoop;
2145 }
2146 
2147 
2148 void Foam::cellCuts::setFromCellLoops
2149 (
2150  const labelList& cellLabels,
2151  const labelListList& cellLoops,
2152  const List<scalarField>& cellLoopWeights
2153 )
2154 {
2155  // 'Uncut' edges/vertices that are not used in loops.
2156  pointIsCut_ = false;
2157 
2158  edgeIsCut_ = false;
2159 
2160  forAll(cellLabels, cellLabelI)
2161  {
2162  label celli = cellLabels[cellLabelI];
2163 
2164  const labelList& loop = cellLoops[cellLabelI];
2165 
2166  if (loop.size())
2167  {
2168  const scalarField& loopWeights = cellLoopWeights[cellLabelI];
2169 
2170  if (setFromCellLoop(celli, loop, loopWeights))
2171  {
2172  // Valid loop. Call above will have upated all already.
2173  }
2174  else
2175  {
2176  // Clear cellLoops
2177  cellLoops_[celli].setSize(0);
2178  }
2179  }
2180  }
2181 }
2182 
2183 
2184 void Foam::cellCuts::setFromCellCutter
2185 (
2186  const cellLooper& cellCutter,
2187  const List<refineCell>& refCells
2188 )
2189 {
2190  // 'Uncut' edges/vertices that are not used in loops.
2191  pointIsCut_ = false;
2192 
2193  edgeIsCut_ = false;
2194 
2195  // storage for loop of cuts (cut vertices and/or cut edges)
2196  labelList cellLoop;
2197  scalarField cellLoopWeights;
2198 
2199  // For debugging purposes
2200  DynamicList<label> invalidCutCells(2);
2201  DynamicList<labelList> invalidCutLoops(2);
2202  DynamicList<scalarField> invalidCutLoopWeights(2);
2203 
2204 
2205  forAll(refCells, refCelli)
2206  {
2207  const refineCell& refCell = refCells[refCelli];
2208 
2209  label celli = refCell.cellNo();
2210 
2211  const vector& refDir = refCell.direction();
2212 
2213 
2214  // Cut cell. Determines cellLoop and cellLoopWeights
2215  bool goodCut =
2216  cellCutter.cut
2217  (
2218  refDir,
2219  celli,
2220 
2221  pointIsCut_,
2222  edgeIsCut_,
2223  edgeWeight_,
2224 
2225  cellLoop,
2226  cellLoopWeights
2227  );
2228 
2229  // Check whether edge refinement is on a per face basis compatible with
2230  // current pattern.
2231  if (goodCut)
2232  {
2233  if (setFromCellLoop(celli, cellLoop, cellLoopWeights))
2234  {
2235  // Valid loop. Will have updated all info already.
2236  }
2237  else
2238  {
2239  cellLoops_[celli].setSize(0);
2240 
2242  << "Found loop on cell " << celli
2243  << " that resulted in an unexpected bad cut." << nl
2244  << " Suggestions:" << nl
2245  << " - Turn on the debug switch for 'cellCuts' to get"
2246  << " geometry files that identify this cell." << nl
2247  << " - Also keep in mind to check the defined"
2248  << " reference directions, as these are most likely the"
2249  << " origin of the problem."
2250  << nl << endl;
2251 
2252  // Discarded by validLoop
2253  if (debug)
2254  {
2255  invalidCutCells.append(celli);
2256  invalidCutLoops.append(cellLoop);
2257  invalidCutLoopWeights.append(cellLoopWeights);
2258  }
2259  }
2260  }
2261  else
2262  {
2263  // Clear cellLoops
2264  cellLoops_[celli].setSize(0);
2265  }
2266  }
2267 
2268  if (debug && invalidCutCells.size())
2269  {
2270  invalidCutCells.shrink();
2271  invalidCutLoops.shrink();
2272  invalidCutLoopWeights.shrink();
2273 
2274  fileName cutsFile("invalidLoopCells.obj");
2275 
2276  Pout<< "cellCuts : writing inValidLoops cells to " << cutsFile << endl;
2277 
2278  OFstream cutsStream(cutsFile);
2279 
2281  (
2282  cutsStream,
2283  mesh().cells(),
2284  mesh().faces(),
2285  mesh().points(),
2286  invalidCutCells
2287  );
2288 
2289  fileName loopsFile("invalidLoops.obj");
2290 
2291  Pout<< "cellCuts : writing inValidLoops loops to " << loopsFile << endl;
2292 
2293  OFstream loopsStream(loopsFile);
2294 
2295  label vertI = 0;
2296 
2297  forAll(invalidCutLoops, i)
2298  {
2299  writeOBJ
2300  (
2301  loopsStream,
2302  loopPoints(invalidCutLoops[i], invalidCutLoopWeights[i]),
2303  vertI
2304  );
2305  }
2306  }
2307 }
2308 
2309 
2310 void Foam::cellCuts::setFromCellCutter
2311 (
2312  const cellLooper& cellCutter,
2313  const labelList& cellLabels,
2314  const List<plane>& cellCutPlanes
2315 )
2316 {
2317  // 'Uncut' edges/vertices that are not used in loops.
2318  pointIsCut_ = false;
2319 
2320  edgeIsCut_ = false;
2321 
2322  // storage for loop of cuts (cut vertices and/or cut edges)
2323  labelList cellLoop;
2324  scalarField cellLoopWeights;
2325 
2326  // For debugging purposes
2327  DynamicList<label> invalidCutCells(2);
2328  DynamicList<labelList> invalidCutLoops(2);
2329  DynamicList<scalarField> invalidCutLoopWeights(2);
2330 
2331 
2332  forAll(cellLabels, i)
2333  {
2334  label celli = cellLabels[i];
2335 
2336  // Cut cell. Determines cellLoop and cellLoopWeights
2337  bool goodCut =
2338  cellCutter.cut
2339  (
2340  cellCutPlanes[i],
2341  celli,
2342 
2343  pointIsCut_,
2344  edgeIsCut_,
2345  edgeWeight_,
2346 
2347  cellLoop,
2348  cellLoopWeights
2349  );
2350 
2351  // Check whether edge refinement is on a per face basis compatible with
2352  // current pattern.
2353  if (goodCut)
2354  {
2355  if (setFromCellLoop(celli, cellLoop, cellLoopWeights))
2356  {
2357  // Valid loop. Will have updated all info already.
2358  }
2359  else
2360  {
2361  cellLoops_[celli].setSize(0);
2362 
2363  // Discarded by validLoop
2364  if (debug)
2365  {
2366  invalidCutCells.append(celli);
2367  invalidCutLoops.append(cellLoop);
2368  invalidCutLoopWeights.append(cellLoopWeights);
2369  }
2370  }
2371  }
2372  else
2373  {
2374  // Clear cellLoops
2375  cellLoops_[celli].setSize(0);
2376  }
2377  }
2378 
2379  if (debug && invalidCutCells.size())
2380  {
2381  invalidCutCells.shrink();
2382  invalidCutLoops.shrink();
2383  invalidCutLoopWeights.shrink();
2384 
2385  fileName cutsFile("invalidLoopCells.obj");
2386 
2387  Pout<< "cellCuts : writing inValidLoops cells to " << cutsFile << endl;
2388 
2389  OFstream cutsStream(cutsFile);
2390 
2392  (
2393  cutsStream,
2394  mesh().cells(),
2395  mesh().faces(),
2396  mesh().points(),
2397  invalidCutCells
2398  );
2399 
2400  fileName loopsFile("invalidLoops.obj");
2401 
2402  Pout<< "cellCuts : writing inValidLoops loops to " << loopsFile << endl;
2403 
2404  OFstream loopsStream(loopsFile);
2405 
2406  label vertI = 0;
2407 
2408  forAll(invalidCutLoops, i)
2409  {
2410  writeOBJ
2411  (
2412  loopsStream,
2413  loopPoints(invalidCutLoops[i], invalidCutLoopWeights[i]),
2414  vertI
2415  );
2416  }
2417  }
2418 }
2419 
2420 
2421 void Foam::cellCuts::orientPlanesAndLoops()
2422 {
2423  // Determine anchorPoints if not yet done by validLoop.
2424  forAll(cellLoops_, celli)
2425  {
2426  const labelList& loop = cellLoops_[celli];
2427 
2428  if (loop.size() && cellAnchorPoints_[celli].empty())
2429  {
2430  // Leave anchor points empty if illegal loop.
2431  calcAnchors
2432  (
2433  celli,
2434  loop,
2435  loopPoints(celli),
2436  cellAnchorPoints_[celli]
2437  );
2438  }
2439  }
2440 
2441  if (debug & 2)
2442  {
2443  Pout<< "cellAnchorPoints:" << endl;
2444  }
2445  forAll(cellAnchorPoints_, celli)
2446  {
2447  if (cellLoops_[celli].size())
2448  {
2449  if (cellAnchorPoints_[celli].empty())
2450  {
2452  << "No anchor points for cut cell " << celli << endl
2453  << "cellLoop:" << cellLoops_[celli] << abort(FatalError);
2454  }
2455 
2456  if (debug & 2)
2457  {
2458  Pout<< " cell:" << celli << " anchored at "
2459  << cellAnchorPoints_[celli] << endl;
2460  }
2461  }
2462  }
2463 
2464  // Calculate number of valid cellLoops
2465  nLoops_ = 0;
2466 
2467  forAll(cellLoops_, celli)
2468  {
2469  if (cellLoops_[celli].size())
2470  {
2471  nLoops_++;
2472  }
2473  }
2474 }
2475 
2476 
2477 void Foam::cellCuts::calcLoopsAndAddressing(const labelList& cutCells)
2478 {
2479  // Sanity check on weights
2480  forAll(edgeIsCut_, edgeI)
2481  {
2482  if (edgeIsCut_[edgeI])
2483  {
2484  scalar weight = edgeWeight_[edgeI];
2485 
2486  if (weight < 0 || weight > 1)
2487  {
2489  << "Weight out of range [0,1]. Edge " << edgeI
2490  << " verts:" << mesh().edges()[edgeI]
2491  << " weight:" << weight << abort(FatalError);
2492  }
2493  }
2494  else
2495  {
2496  // Weight not used. Set to illegal value to make any use fall over.
2497  edgeWeight_[edgeI] = -GREAT;
2498  }
2499  }
2500 
2501 
2502  // Calculate faces that split cells in two
2503  calcCellLoops(cutCells);
2504 
2505  if (debug & 2)
2506  {
2507  Pout<< "-- cellLoops --" << endl;
2508  forAll(cellLoops_, celli)
2509  {
2510  const labelList& loop = cellLoops_[celli];
2511 
2512  if (loop.size())
2513  {
2514  Pout<< "cell:" << celli << " ";
2515  writeCuts(Pout, loop, loopWeights(loop));
2516  Pout<< endl;
2517  }
2518  }
2519  }
2520 
2521  // Redo basic cut information (pointIsCut, edgeIsCut, faceSplitCut)
2522  // using cellLoop only.
2523  setFromCellLoops();
2524 }
2525 
2526 
2527 void Foam::cellCuts::check() const
2528 {
2529  // Check weights for unsnapped values
2530  forAll(edgeIsCut_, edgeI)
2531  {
2532  if (edgeIsCut_[edgeI])
2533  {
2534  if
2535  (
2536  edgeWeight_[edgeI] < geomCellLooper::snapTol()
2537  || edgeWeight_[edgeI] > (1 - geomCellLooper::snapTol())
2538  )
2539  {
2540  // Should have been snapped.
2542  << "edge:" << edgeI << " vertices:"
2543  << mesh().edges()[edgeI]
2544  << " weight:" << edgeWeight_[edgeI] << " should have been"
2545  << " snapped to one of its endpoints"
2546  << endl;
2547  }
2548  }
2549  else
2550  {
2551  if (edgeWeight_[edgeI] > - 1)
2552  {
2554  << "edge:" << edgeI << " vertices:"
2555  << mesh().edges()[edgeI]
2556  << " weight:" << edgeWeight_[edgeI] << " is not cut but"
2557  << " its weight is not marked invalid"
2558  << abort(FatalError);
2559  }
2560  }
2561  }
2562 
2563  // Check that all elements of cellloop are registered
2564  forAll(cellLoops_, celli)
2565  {
2566  const labelList& loop = cellLoops_[celli];
2567 
2568  forAll(loop, i)
2569  {
2570  label cut = loop[i];
2571 
2572  if
2573  (
2574  (isEdge(cut) && !edgeIsCut_[getEdge(cut)])
2575  || (!isEdge(cut) && !pointIsCut_[getVertex(cut)])
2576  )
2577  {
2578  labelList cuts(1, cut);
2579  writeCuts(Pout, cuts, loopWeights(cuts));
2580 
2582  << "cell:" << celli << " loop:"
2583  << loop
2584  << " cut:" << cut << " is not marked as cut"
2585  << abort(FatalError);
2586  }
2587  }
2588  }
2589 
2590  // Check that no elements of cell loop are anchor point.
2591  forAll(cellLoops_, celli)
2592  {
2593  const labelList& anchors = cellAnchorPoints_[celli];
2594 
2595  const labelList& loop = cellLoops_[celli];
2596 
2597  if (loop.size() && anchors.empty())
2598  {
2600  << "cell:" << celli << " loop:" << loop
2601  << " has no anchor points"
2602  << abort(FatalError);
2603  }
2604 
2605 
2606  forAll(loop, i)
2607  {
2608  label cut = loop[i];
2609 
2610  if
2611  (
2612  !isEdge(cut)
2613  && findIndex(anchors, getVertex(cut)) != -1
2614  )
2615  {
2617  << "cell:" << celli << " loop:" << loop
2618  << " anchor points:" << anchors
2619  << " anchor:" << getVertex(cut) << " is part of loop"
2620  << abort(FatalError);
2621  }
2622  }
2623  }
2624 
2625 
2626  // Check that cut faces have a neighbour that is cut.
2627  forAllConstIter(Map<edge>, faceSplitCut_, iter)
2628  {
2629  label facei = iter.key();
2630 
2631  if (mesh().isInternalFace(facei))
2632  {
2633  label own = mesh().faceOwner()[facei];
2634  label nei = mesh().faceNeighbour()[facei];
2635 
2636  if (cellLoops_[own].empty() && cellLoops_[nei].empty())
2637  {
2639  << "Internal face:" << facei << " cut by " << iter()
2640  << " has owner:" << own
2641  << " and neighbour:" << nei
2642  << " that are both uncut"
2643  << abort(FatalError);
2644  }
2645  }
2646  else
2647  {
2648  label own = mesh().faceOwner()[facei];
2649 
2650  if (cellLoops_[own].empty())
2651  {
2653  << "Boundary face:" << facei << " cut by " << iter()
2654  << " has owner:" << own
2655  << " that is uncut"
2656  << abort(FatalError);
2657  }
2658  }
2659  }
2660 }
2661 
2662 
2663 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
2664 
2665 Foam::cellCuts::cellCuts
2667  const polyMesh& mesh,
2668  const labelList& cutCells,
2669  const labelList& meshVerts,
2670  const labelList& meshEdges,
2671  const scalarField& meshEdgeWeights
2672 )
2673 :
2674  edgeVertex(mesh),
2675  pointIsCut_(expand(mesh.nPoints(), meshVerts)),
2676  edgeIsCut_(expand(mesh.nEdges(), meshEdges)),
2677  edgeWeight_(expand(mesh.nEdges(), meshEdges, meshEdgeWeights)),
2678  faceCutsPtr_(NULL),
2679  faceSplitCut_(cutCells.size()),
2680  cellLoops_(mesh.nCells()),
2681  nLoops_(-1),
2682  cellAnchorPoints_(mesh.nCells())
2683 {
2684  if (debug)
2685  {
2686  Pout<< "cellCuts : constructor from cut verts and edges" << endl;
2687  }
2688 
2689  calcLoopsAndAddressing(cutCells);
2690 
2691  // Calculate planes and flip cellLoops if necessary
2692  orientPlanesAndLoops();
2693 
2694  if (debug)
2695  {
2696  check();
2697  }
2698 
2699  clearOut();
2700 
2701  if (debug)
2702  {
2703  Pout<< "cellCuts : leaving constructor from cut verts and edges"
2704  << endl;
2705  }
2706 }
2707 
2708 
2709 Foam::cellCuts::cellCuts
2711  const polyMesh& mesh,
2712  const labelList& meshVerts,
2713  const labelList& meshEdges,
2714  const scalarField& meshEdgeWeights
2715 )
2716 :
2717  edgeVertex(mesh),
2718  pointIsCut_(expand(mesh.nPoints(), meshVerts)),
2719  edgeIsCut_(expand(mesh.nEdges(), meshEdges)),
2720  edgeWeight_(expand(mesh.nEdges(), meshEdges, meshEdgeWeights)),
2721  faceCutsPtr_(NULL),
2722  faceSplitCut_(mesh.nFaces()/10 + 1),
2723  cellLoops_(mesh.nCells()),
2724  nLoops_(-1),
2725  cellAnchorPoints_(mesh.nCells())
2726 {
2727  // Construct from pattern of cuts. Finds out itself which cells are cut.
2728  // (can go wrong if e.g. all neighbours of cell are refined)
2729 
2730  if (debug)
2731  {
2732  Pout<< "cellCuts : constructor from cellLoops" << endl;
2733  }
2734 
2735  calcLoopsAndAddressing(identity(mesh.nCells()));
2736 
2737  // Calculate planes and flip cellLoops if necessary
2738  orientPlanesAndLoops();
2739 
2740  if (debug)
2741  {
2742  check();
2743  }
2744 
2745  clearOut();
2746 
2747  if (debug)
2748  {
2749  Pout<< "cellCuts : leaving constructor from cellLoops" << endl;
2750  }
2751 }
2752 
2753 
2754 Foam::cellCuts::cellCuts
2756  const polyMesh& mesh,
2757  const labelList& cellLabels,
2758  const labelListList& cellLoops,
2759  const List<scalarField>& cellEdgeWeights
2760 )
2761 :
2762  edgeVertex(mesh),
2763  pointIsCut_(mesh.nPoints(), false),
2764  edgeIsCut_(mesh.nEdges(), false),
2765  edgeWeight_(mesh.nEdges(), -GREAT),
2766  faceCutsPtr_(NULL),
2767  faceSplitCut_(cellLabels.size()),
2768  cellLoops_(mesh.nCells()),
2769  nLoops_(-1),
2770  cellAnchorPoints_(mesh.nCells())
2771 {
2772  if (debug)
2773  {
2774  Pout<< "cellCuts : constructor from cellLoops" << endl;
2775  }
2776 
2777  // Update pointIsCut, edgeIsCut, faceSplitCut from cell loops.
2778  // Makes sure cuts are consistent
2779  setFromCellLoops(cellLabels, cellLoops, cellEdgeWeights);
2780 
2781  // Calculate planes and flip cellLoops if necessary
2782  orientPlanesAndLoops();
2783 
2784  if (debug)
2785  {
2786  check();
2787  }
2788 
2789  clearOut();
2790 
2791  if (debug)
2792  {
2793  Pout<< "cellCuts : leaving constructor from cellLoops" << endl;
2794  }
2795 }
2796 
2797 
2798 Foam::cellCuts::cellCuts
2800  const polyMesh& mesh,
2801  const cellLooper& cellCutter,
2802  const List<refineCell>& refCells
2803 )
2804 :
2805  edgeVertex(mesh),
2806  pointIsCut_(mesh.nPoints(), false),
2807  edgeIsCut_(mesh.nEdges(), false),
2808  edgeWeight_(mesh.nEdges(), -GREAT),
2809  faceCutsPtr_(NULL),
2810  faceSplitCut_(refCells.size()),
2811  cellLoops_(mesh.nCells()),
2812  nLoops_(-1),
2813  cellAnchorPoints_(mesh.nCells())
2814 {
2815  if (debug)
2816  {
2817  Pout<< "cellCuts : constructor from cellCutter" << endl;
2818  }
2819 
2820  // Update pointIsCut, edgeIsCut, faceSplitCut from cell loops.
2821  // Makes sure cuts are consistent
2822  setFromCellCutter(cellCutter, refCells);
2823 
2824  // Calculate planes and flip cellLoops if necessary
2825  orientPlanesAndLoops();
2826 
2827  if (debug)
2828  {
2829  check();
2830  }
2831 
2832  clearOut();
2833 
2834  if (debug)
2835  {
2836  Pout<< "cellCuts : leaving constructor from cellCutter" << endl;
2837  }
2838 }
2839 
2840 
2841 Foam::cellCuts::cellCuts
2843  const polyMesh& mesh,
2844  const cellLooper& cellCutter,
2845  const labelList& cellLabels,
2846  const List<plane>& cutPlanes
2847 )
2848 :
2849  edgeVertex(mesh),
2850  pointIsCut_(mesh.nPoints(), false),
2851  edgeIsCut_(mesh.nEdges(), false),
2852  edgeWeight_(mesh.nEdges(), -GREAT),
2853  faceCutsPtr_(NULL),
2854  faceSplitCut_(cellLabels.size()),
2855  cellLoops_(mesh.nCells()),
2856  nLoops_(-1),
2857  cellAnchorPoints_(mesh.nCells())
2858 {
2859  if (debug)
2860  {
2861  Pout<< "cellCuts : constructor from cellCutter with prescribed plane"
2862  << endl;
2863  }
2864 
2865  // Update pointIsCut, edgeIsCut, faceSplitCut from cell loops.
2866  // Makes sure cuts are consistent
2867  setFromCellCutter(cellCutter, cellLabels, cutPlanes);
2868 
2869  // Calculate planes and flip cellLoops if necessary
2870  orientPlanesAndLoops();
2871 
2872  if (debug)
2873  {
2874  check();
2875  }
2876 
2877  clearOut();
2878 
2879  if (debug)
2880  {
2881  Pout<< "cellCuts : leaving constructor from cellCutter with prescribed"
2882  << " plane" << endl;
2883  }
2884 }
2885 
2886 
2887 Foam::cellCuts::cellCuts
2889  const polyMesh& mesh,
2890  const boolList& pointIsCut,
2891  const boolList& edgeIsCut,
2892  const scalarField& edgeWeight,
2893  const Map<edge>& faceSplitCut,
2894  const labelListList& cellLoops,
2895  const label nLoops,
2896  const labelListList& cellAnchorPoints
2897 )
2898 :
2899  edgeVertex(mesh),
2900  pointIsCut_(pointIsCut),
2901  edgeIsCut_(edgeIsCut),
2902  edgeWeight_(edgeWeight),
2903  faceCutsPtr_(NULL),
2904  faceSplitCut_(faceSplitCut),
2905  cellLoops_(cellLoops),
2906  nLoops_(nLoops),
2907  cellAnchorPoints_(cellAnchorPoints)
2908 {
2909  if (debug)
2910  {
2911  Pout<< "cellCuts : constructor from components" << endl;
2912  Pout<< "cellCuts : leaving constructor from components" << endl;
2913  }
2914 }
2915 
2916 
2917 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
2918 
2920 {
2921  clearOut();
2922 }
2923 
2924 
2926 {
2927  deleteDemandDrivenData(faceCutsPtr_);
2928 }
2929 
2930 
2931 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
2932 
2933 Foam::pointField Foam::cellCuts::loopPoints(const label celli) const
2934 {
2935  const labelList& loop = cellLoops_[celli];
2936 
2937  pointField loopPts(loop.size());
2938 
2939  forAll(loop, fp)
2940  {
2941  label cut = loop[fp];
2942 
2943  if (isEdge(cut))
2944  {
2945  label edgeI = getEdge(cut);
2946 
2947  loopPts[fp] = coord(cut, edgeWeight_[edgeI]);
2948  }
2949  else
2950  {
2951  loopPts[fp] = coord(cut, -GREAT);
2952  }
2953  }
2954  return loopPts;
2955 }
2956 
2957 
2958 void Foam::cellCuts::flip(const label celli)
2959 {
2960  labelList& loop = cellLoops_[celli];
2961 
2962  reverse(loop);
2963 
2964  // Reverse anchor point set.
2965  cellAnchorPoints_[celli] =
2966  nonAnchorPoints
2967  (
2968  mesh().cellPoints()[celli],
2969  cellAnchorPoints_[celli],
2970  loop
2971  );
2972 }
2973 
2974 
2976 {
2977  labelList& loop = cellLoops_[celli];
2978 
2979  reverse(loop);
2980 }
2981 
2982 
2983 void Foam::cellCuts::writeOBJ
2985  Ostream& os,
2986  const pointField& loopPts,
2987  label& vertI
2988 ) const
2989 {
2990  label startVertI = vertI;
2991 
2992  forAll(loopPts, fp)
2993  {
2994  const point& pt = loopPts[fp];
2995 
2996  os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
2997 
2998  vertI++;
2999  }
3000 
3001  os << 'f';
3002  forAll(loopPts, fp)
3003  {
3004  os << ' ' << startVertI + fp + 1;
3005  }
3006  os << endl;
3007 }
3008 
3009 
3010 void Foam::cellCuts::writeOBJ(Ostream& os) const
3011 {
3012  label vertI = 0;
3013 
3014  forAll(cellLoops_, celli)
3015  {
3016  writeOBJ(os, loopPoints(celli), vertI);
3017  }
3018 }
3019 
3020 
3021 void Foam::cellCuts::writeCellOBJ(const fileName& dir, const label celli) const
3022 {
3023  const labelList& anchors = cellAnchorPoints_[celli];
3024 
3025  writeOBJ(dir, celli, loopPoints(celli), anchors);
3026 }
3027 
3028 
3029 // ************************************************************************* //
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
point centre(const pointField &) const
Centre point of face.
Definition: face.C:488
const Cmpt & z() const
Definition: VectorI.H:87
#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 class for handling file names.
Definition: fileName.H:69
An STL-conforming const_iterator.
Definition: HashTable.H:470
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const Cmpt & x() const
Definition: VectorI.H:75
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:313
bool edgeOnCell(const primitiveMesh &, const label celli, const label edgeI)
Is edge used by cell.
Definition: meshTools.C:285
Output to file stream.
Definition: OFstream.H:81
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
const labelListList & cellPoints() const
List< face > faceList
Definition: faceListFwd.H:43
void flip(const label celli)
Flip loop for celli. Updates anchor points as well.
Definition: cellCuts.C:2958
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Abstract base class. Concrete implementations know how to cut a cell (i.e. determine a loop around th...
Definition: cellLooper.H:69
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
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
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
vector normal(const pointField &) const
Vector normal; magnitude is equal to area of face.
Definition: face.C:552
const cellList & cells() const
bool faceOnCell(const primitiveMesh &, const label celli, const label facei)
Is face used by cell.
Definition: meshTools.C:308
scalar f1
Definition: createFields.H:28
const labelListList & edgeFaces() const
const vector & direction() const
Definition: refineCell.H:87
label nCells() const
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Various functions to operate on Lists.
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
bool insert(const label &, const T &newElmt)
Insert a new hashedEntry.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
static const labelSphericalTensor labelI(1)
Identity labelTensor.
dynamicFvMesh & mesh
const cellShapeList & cells
const pointField & points
List< edge > edgeList
Definition: edgeList.H:38
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
label cellNo() const
Definition: refineCell.H:82
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void clearOut()
Clear out demand driven storage.
Definition: cellCuts.C:2925
string expand(const string &, const HashTable< string, word, string::hash > &mapping, const char sigil= '$')
Expand occurences of variables according to the mapping.
Definition: stringOps.C:75
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
const Cmpt & y() const
Definition: VectorI.H:81
Combines edge or vertex in single label. Used to specify cuts across cell circumference.
Definition: edgeVertex.H:52
List< label > labelList
A List of labels.
Definition: labelList.H:56
const labelListList & pointEdges() const
static const zero Zero
Definition: zero.H:91
Container with cells to refine. Refinement given as single direction.
Definition: refineCell.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
~cellCuts()
Destructor.
Definition: cellCuts.C:2919
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
label start() const
Return start vertex label.
Definition: edgeI.H:81
void reverse(UList< T > &, const label n)
Definition: UListI.H:322
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:240
static const char nl
Definition: Ostream.H:262
const labelListList & cellEdges() const
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
label nEdges() const
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
A normal distribution model.
label nFaces() const
void setSize(const label)
Reset size of List.
Definition: List.C:295
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
#define WarningInFunction
Report a warning using Foam::Warning.
virtual bool cut(const vector &refDir, const label celli, const boolList &vertIsCut, const boolList &edgeIsCut, const scalarField &edgeWeight, labelList &loop, scalarField &loopWeights) const =0
Create cut along circumference of celli. Gets current mesh cuts.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
label nPoints() const
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
label walkFace(const primitiveMesh &, const label facei, const label startEdgeI, const label startVertI, const label nEdges)
Returns label of edge nEdges away from startEdge (in the direction.
Definition: meshTools.C:587
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
const labelListList & faceEdges() const
void deleteDemandDrivenData(DataPtr &dataPtr)
void writeCellOBJ(const fileName &dir, const label celli) const
debugging:Write edges of cell and loop
Definition: cellCuts.C:3021
static scalar snapTol()
labelList nonAnchorPoints(const labelList &cellPoints, const labelList &anchorPoints, const labelList &loop) const
Invert anchor point selection.
Definition: cellCuts.C:1167
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:365
void getEdgeFaces(const primitiveMesh &, const label celli, const label edgeI, label &face0, label &face1)
Get faces on cell using edgeI. Throws error if no two found.
Definition: meshTools.C:458
static const label labelMin
Definition: label.H:61
label otherFace(const primitiveMesh &, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:536
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
void flipLoopOnly(const label celli)
Flip loop for celli. Does not update anchors. Use with care.
Definition: cellCuts.C:2975