cellClassification.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 "cellClassification.H"
27 #include "triSurfaceSearch.H"
28 #include "indexedOctree.H"
29 #include "treeDataFace.H"
30 #include "meshSearch.H"
31 #include "cellInfo.H"
32 #include "polyMesh.H"
33 #include "FaceCellWave.H"
34 #include "ListOps.H"
35 #include "meshTools.H"
36 #include "cpuTime.H"
37 #include "triSurface.H"
38 #include "globalMeshData.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 defineTypeNameAndDebug(cellClassification, 0);
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 Foam::label Foam::cellClassification::count
51 (
52  const labelList& elems,
53  const label elem
54 )
55 {
56  label cnt = 0;
57 
58  forAll(elems, elemI)
59  {
60  if (elems[elemI] == elem)
61  {
62  cnt++;
63  }
64  }
65  return cnt;
66 
67 }
68 
69 
70 // Mark all faces that are cut by the surface. Two pass:
71 // Pass1: mark all mesh edges that intersect surface (quick since triangle
72 // pierce test).
73 // Pass2: Check for all point neighbours of these faces whether any of their
74 // faces are pierced.
75 Foam::boolList Foam::cellClassification::markFaces
76 (
77  const triSurfaceSearch& search
78 ) const
79 {
80  cpuTime timer;
81 
82  boolList cutFace(mesh_.nFaces(), false);
83 
84  label nCutFaces = 0;
85 
86  // Intersect mesh edges with surface (is fast) and mark all faces that
87  // use edge.
88 
89  forAll(mesh_.edges(), edgeI)
90  {
91  if (debug && (edgeI % 10000 == 0))
92  {
93  Pout<< "Intersecting mesh edge " << edgeI << " with surface"
94  << endl;
95  }
96 
97  const edge& e = mesh_.edges()[edgeI];
98 
99  const point& p0 = mesh_.points()[e.start()];
100  const point& p1 = mesh_.points()[e.end()];
101 
102  pointIndexHit pHit(search.tree().findLineAny(p0, p1));
103 
104  if (pHit.hit())
105  {
106  const labelList& myFaces = mesh_.edgeFaces()[edgeI];
107 
108  forAll(myFaces, myFacei)
109  {
110  label facei = myFaces[myFacei];
111 
112  if (!cutFace[facei])
113  {
114  cutFace[facei] = true;
115 
116  nCutFaces++;
117  }
118  }
119  }
120  }
121 
122  if (debug)
123  {
124  Pout<< "Intersected edges of mesh with surface in = "
125  << timer.cpuTimeIncrement() << " s\n" << endl << endl;
126  }
127 
128  //
129  // Construct octree on faces that have not yet been marked as cut
130  //
131 
132  labelList allFaces(mesh_.nFaces() - nCutFaces);
133 
134  label allFacei = 0;
135 
136  forAll(cutFace, facei)
137  {
138  if (!cutFace[facei])
139  {
140  allFaces[allFacei++] = facei;
141  }
142  }
143 
144  if (debug)
145  {
146  Pout<< "Testing " << allFacei << " faces for piercing by surface"
147  << endl;
148  }
149 
150  treeBoundBox allBb(mesh_.points());
151 
152  // Extend domain slightly (also makes it 3D if was 2D)
153  scalar tol = 1e-6 * allBb.avgDim();
154 
155  point& bbMin = allBb.min();
156  bbMin.x() -= tol;
157  bbMin.y() -= tol;
158  bbMin.z() -= tol;
159 
160  point& bbMax = allBb.max();
161  bbMax.x() += 2*tol;
162  bbMax.y() += 2*tol;
163  bbMax.z() += 2*tol;
164 
165  indexedOctree<treeDataFace> faceTree
166  (
167  treeDataFace(false, mesh_, allFaces),
168  allBb, // overall search domain
169  8, // maxLevel
170  10, // leafsize
171  3.0 // duplicity
172  );
173 
174  const triSurface& surf = search.surface();
175  const edgeList& edges = surf.edges();
176  const pointField& localPoints = surf.localPoints();
177 
178  label nAddFaces = 0;
179 
180  forAll(edges, edgeI)
181  {
182  if (debug && (edgeI % 10000 == 0))
183  {
184  Pout<< "Intersecting surface edge " << edgeI
185  << " with mesh faces" << endl;
186  }
187  const edge& e = edges[edgeI];
188 
189  const point& start = localPoints[e.start()];
190  const point& end = localPoints[e.end()];
191 
192  vector edgeNormal(end - start);
193  const scalar edgeMag = mag(edgeNormal);
194  const vector smallVec = 1e-9*edgeNormal;
195 
196  edgeNormal /= edgeMag+vSmall;
197 
198  // Current start of pierce test
199  point pt = start;
200 
201  while (true)
202  {
203  pointIndexHit pHit(faceTree.findLine(pt, end));
204 
205  if (!pHit.hit())
206  {
207  break;
208  }
209  else
210  {
211  label facei = faceTree.shapes().faceLabels()[pHit.index()];
212 
213  if (!cutFace[facei])
214  {
215  cutFace[facei] = true;
216 
217  nAddFaces++;
218  }
219 
220  // Restart from previous endpoint
221  pt = pHit.hitPoint() + smallVec;
222 
223  if (((pt-start) & edgeNormal) >= edgeMag)
224  {
225  break;
226  }
227  }
228  }
229  }
230 
231  if (debug)
232  {
233  Pout<< "Detected an additional " << nAddFaces << " faces cut"
234  << endl;
235 
236  Pout<< "Intersected edges of surface with mesh faces in = "
237  << timer.cpuTimeIncrement() << " s\n" << endl << endl;
238  }
239 
240  return cutFace;
241 }
242 
243 
244 // Determine faces cut by surface and use to divide cells into types. See
245 // cellInfo. All cells reachable from outsidePts are considered to be of type
246 // 'outside'
247 void Foam::cellClassification::markCells
248 (
249  const meshSearch& queryMesh,
250  const boolList& piercedFace,
251  const pointField& outsidePts
252 )
253 {
254  // Use meshwave to partition mesh, starting from outsidePt
255 
256 
257  // Construct null; sets type to NOTSET
258  List<cellInfo> cellInfoList(mesh_.nCells());
259 
260  // Mark cut cells first
261  forAll(piercedFace, facei)
262  {
263  if (piercedFace[facei])
264  {
265  cellInfoList[mesh_.faceOwner()[facei]] =
266  cellInfo(cellClassification::CUT);
267 
268  if (mesh_.isInternalFace(facei))
269  {
270  cellInfoList[mesh_.faceNeighbour()[facei]] =
271  cellInfo(cellClassification::CUT);
272  }
273  }
274  }
275 
276  //
277  // Mark cells containing outside points as being outside
278  //
279 
280  // Coarse guess number of faces
281  labelHashSet outsideFacesMap(outsidePts.size() * 6 * 2);
282 
283  forAll(outsidePts, outsidePtI)
284  {
285  // Use linear search for points.
286  label celli = queryMesh.findCell(outsidePts[outsidePtI], -1, false);
287 
288  if (returnReduce(celli, maxOp<label>()) == -1)
289  {
291  << "outsidePoint " << outsidePts[outsidePtI]
292  << " is not inside any cell"
293  << nl << "It might be on a face or outside the geometry"
294  << exit(FatalError);
295  }
296 
297  if (celli >= 0)
298  {
299  cellInfoList[celli] = cellInfo(cellClassification::OUTSIDE);
300 
301  // Mark faces of celli
302  const labelList& myFaces = mesh_.cells()[celli];
303  forAll(myFaces, myFacei)
304  {
305  outsideFacesMap.insert(myFaces[myFacei]);
306  }
307  }
308  }
309 
310 
311  //
312  // Mark faces to start wave from
313  //
314 
315  labelList changedFaces(outsideFacesMap.toc());
316 
317  List<cellInfo> changedFacesInfo
318  (
319  changedFaces.size(),
321  );
322 
323  List<cellInfo> faceInfoList(mesh().nFaces());
324  FaceCellWave<cellInfo> cellInfoCalc
325  (
326  mesh_,
327  changedFaces, // Labels of changed faces
328  changedFacesInfo, // Information on changed faces
329  faceInfoList,
330  cellInfoList, // Information on all cells
331  mesh_.globalData().nTotalCells() + 1 // max iterations
332  );
333 
334  forAll(cellInfoList, celli)
335  {
336  label t = cellInfoList[celli].type();
337 
339  {
341  }
342  operator[](celli) = t;
343  }
344 }
345 
346 
347 void Foam::cellClassification::classifyPoints
348 (
349  const label meshType,
350  const labelList& cellType,
351  List<pointStatus>& pointSide
352 ) const
353 {
354  pointSide.setSize(mesh_.nPoints());
355 
356  forAll(mesh_.pointCells(), pointi)
357  {
358  const labelList& pCells = mesh_.pointCells()[pointi];
359 
360  pointSide[pointi] = UNSET;
361 
362  forAll(pCells, i)
363  {
364  label type = cellType[pCells[i]];
365 
366  if (type == meshType)
367  {
368  if (pointSide[pointi] == UNSET)
369  {
370  pointSide[pointi] = MESH;
371  }
372  else if (pointSide[pointi] == NONMESH)
373  {
374  pointSide[pointi] = MIXED;
375 
376  break;
377  }
378  }
379  else
380  {
381  if (pointSide[pointi] == UNSET)
382  {
383  pointSide[pointi] = NONMESH;
384  }
385  else if (pointSide[pointi] == MESH)
386  {
387  pointSide[pointi] = MIXED;
388 
389  break;
390  }
391  }
392  }
393  }
394 }
395 
396 
397 bool Foam::cellClassification::usesMixedPointsOnly
398 (
399  const List<pointStatus>& pointSide,
400  const label celli
401 ) const
402 {
403  const faceList& faces = mesh_.faces();
404 
405  const cell& cFaces = mesh_.cells()[celli];
406 
407  forAll(cFaces, cFacei)
408  {
409  const face& f = faces[cFaces[cFacei]];
410 
411  forAll(f, fp)
412  {
413  if (pointSide[f[fp]] != MIXED)
414  {
415  return false;
416  }
417  }
418  }
419 
420  // All points are mixed.
421  return true;
422 }
423 
424 
425 void Foam::cellClassification::getMeshOutside
426 (
427  const label meshType,
428  faceList& outsideFaces,
429  labelList& outsideOwner
430 ) const
431 {
432  const labelList& own = mesh_.faceOwner();
433  const labelList& nbr = mesh_.faceNeighbour();
434  const faceList& faces = mesh_.faces();
435 
436  outsideFaces.setSize(mesh_.nFaces());
437  outsideOwner.setSize(mesh_.nFaces());
438  label outsideI = 0;
439 
440  // Get faces on interface between meshType and non-meshType
441 
442  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
443  {
444  label ownType = operator[](own[facei]);
445  label nbrType = operator[](nbr[facei]);
446 
447  if (ownType == meshType && nbrType != meshType)
448  {
449  outsideFaces[outsideI] = faces[facei];
450  outsideOwner[outsideI] = own[facei]; // meshType cell
451  outsideI++;
452  }
453  else if (ownType != meshType && nbrType == meshType)
454  {
455  outsideFaces[outsideI] = faces[facei];
456  outsideOwner[outsideI] = nbr[facei]; // meshType cell
457  outsideI++;
458  }
459  }
460 
461  // Get faces on outside of real mesh with cells of meshType.
462 
463  for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
464  {
465  if (operator[](own[facei]) == meshType)
466  {
467  outsideFaces[outsideI] = faces[facei];
468  outsideOwner[outsideI] = own[facei]; // meshType cell
469  outsideI++;
470  }
471  }
472  outsideFaces.setSize(outsideI);
473  outsideOwner.setSize(outsideI);
474 }
475 
476 
477 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
478 
479 // Construct from mesh and surface and point(s) on outside
481 (
482  const polyMesh& mesh,
483  const meshSearch& meshQuery,
484  const triSurfaceSearch& surfQuery,
485  const pointField& outsidePoints
486 )
487 :
489  mesh_(mesh)
490 {
491  markCells
492  (
493  meshQuery,
494  markFaces(surfQuery),
495  outsidePoints
496  );
497 }
498 
499 
500 // Construct from mesh and meshType.
502 (
503  const polyMesh& mesh,
504  const labelList& cellType
505 )
506 :
507  labelList(cellType),
508  mesh_(mesh)
509 {
510  if (mesh_.nCells() != size())
511  {
513  << "Number of elements of cellType argument is not equal to the"
514  << " number of cells" << abort(FatalError);
515  }
516 }
517 
518 
519 // Copy constructor
521 :
522  labelList(cType),
523  mesh_(cType.mesh())
524 {}
525 
526 
527 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
528 
529 
530 // Makes cutCells further than nLayers away from meshType to
531 // fillType. Returns number of cells changed.
533 (
534  const label nLayers,
535  const label meshType,
536  const label fillType
537 )
538 {
539  // Temporary cell type for growing.
540  labelList newCellType(*this);
541 
542 // // Split types into outside and rest
543 // forAll(*this, celli)
544 // {
545 // label type = operator[](celli);
546 //
547 // if (type == meshType)
548 // {
549 // newCellType[celli] = type;
550 // }
551 // else
552 // {
553 // newCellType[celli] = fillType;
554 // }
555 // }
556 
557  newCellType = *this;
558 
559 
560  // Do point-cell-point walk on newCellType out from meshType cells
561  for (label iter = 0; iter < nLayers; iter++)
562  {
563  // Get status of points: visible from meshType (pointSide == MESH)
564  // or fillType/cutCells cells (pointSide == NONMESH) or
565  // both (pointSide == MIXED)
566  List<pointStatus> pointSide(mesh_.nPoints());
567  classifyPoints(meshType, newCellType, pointSide);
568 
569  // Grow layer of outside cells
570  forAll(pointSide, pointi)
571  {
572  if (pointSide[pointi] == MIXED)
573  {
574  // Make cut
575  const labelList& pCells = mesh_.pointCells()[pointi];
576 
577  forAll(pCells, i)
578  {
579  label type = newCellType[pCells[i]];
580 
581  if (type == cellClassification::CUT)
582  {
583  // Found cut cell sharing point with
584  // mesh type cell. Grow.
585  newCellType[pCells[i]] = meshType;
586  }
587  }
588  }
589  }
590  }
591 
592  // Merge newCellType into *this:
593  // - leave meshType cells intact
594  // - leave nonMesh cells intact
595  // - make cutcells fillType if they were not marked in newCellType
596 
597  label nChanged = 0;
598 
599  forAll(newCellType, celli)
600  {
601  if (operator[](celli) == cellClassification::CUT)
602  {
603  if (newCellType[celli] != meshType)
604  {
605  // Cell was cutCell but further than nLayers away from
606  // meshType. Convert to fillType.
607  operator[](celli) = fillType;
608  nChanged++;
609  }
610  }
611  }
612 
613  return nChanged;
614 }
615 
616 
617 // Grow surface by pointNeighbours of type meshType
619 (
620  const label meshType,
621  const label fillType
622 )
623 {
624  boolList hasMeshType(mesh_.nPoints(), false);
625 
626  // Mark points used by meshType cells
627 
628  forAll(mesh_.pointCells(), pointi)
629  {
630  const labelList& myCells = mesh_.pointCells()[pointi];
631 
632  // Check if one of cells has meshType
633  forAll(myCells, myCelli)
634  {
635  label type = operator[](myCells[myCelli]);
636 
637  if (type == meshType)
638  {
639  hasMeshType[pointi] = true;
640 
641  break;
642  }
643  }
644  }
645 
646  // Change neighbours of marked points
647 
648  label nChanged = 0;
649 
650  forAll(hasMeshType, pointi)
651  {
652  if (hasMeshType[pointi])
653  {
654  const labelList& myCells = mesh_.pointCells()[pointi];
655 
656  forAll(myCells, myCelli)
657  {
658  if (operator[](myCells[myCelli]) != meshType)
659  {
660  operator[](myCells[myCelli]) = fillType;
661 
662  nChanged++;
663  }
664  }
665  }
666  }
667  return nChanged;
668 }
669 
670 
671 // Check all points used by cells of meshType for use of at least one point
672 // which is not on the outside. If all points are on the outside of the mesh
673 // this would probably result in a flattened cell when projecting the cell
674 // onto the surface.
676 (
677  const label meshType,
678  const label fillType,
679  const label maxIter
680 )
681 {
682  label nTotChanged = 0;
683 
684  for (label iter = 0; iter < maxIter; iter++)
685  {
686  label nChanged = 0;
687 
688  // Get status of points: visible from meshType or non-meshType cells
689  List<pointStatus> pointSide(mesh_.nPoints());
690  classifyPoints(meshType, *this, pointSide);
691 
692  // Check all cells using mixed point type for whether they use mixed
693  // points only. Note: could probably speed this up by counting number
694  // of mixed verts per face and mixed faces per cell or something?
695  forAll(pointSide, pointi)
696  {
697  if (pointSide[pointi] == MIXED)
698  {
699  const labelList& pCells = mesh_.pointCells()[pointi];
700 
701  forAll(pCells, i)
702  {
703  label celli = pCells[i];
704 
705  if (operator[](celli) == meshType)
706  {
707  if (usesMixedPointsOnly(pointSide, celli))
708  {
709  operator[](celli) = fillType;
710 
711  nChanged++;
712  }
713  }
714  }
715  }
716  }
717  nTotChanged += nChanged;
718 
719  Pout<< "removeHangingCells : changed " << nChanged
720  << " hanging cells" << endl;
721 
722  if (nChanged == 0)
723  {
724  break;
725  }
726  }
727 
728  return nTotChanged;
729 }
730 
731 
733 (
734  const label meshType,
735  const label fillType,
736  const label maxIter
737 )
738 {
739  label nTotChanged = 0;
740 
741  for (label iter = 0; iter < maxIter; iter++)
742  {
743  // Get interface between meshType cells and non-meshType cells as a list
744  // of faces and for each face the cell which is the meshType.
745  faceList outsideFaces;
746  labelList outsideOwner;
747 
748  getMeshOutside(meshType, outsideFaces, outsideOwner);
749 
750  // Build primitivePatch out of it and check it for problems.
751  primitiveFacePatch fp(outsideFaces, mesh_.points());
752 
753  const labelListList& edgeFaces = fp.edgeFaces();
754 
755  label nChanged = 0;
756 
757  // Check all edgeFaces for non-manifoldness
758 
759  forAll(edgeFaces, edgeI)
760  {
761  const labelList& eFaces = edgeFaces[edgeI];
762 
763  if (eFaces.size() > 2)
764  {
765  // patch connected through pinched edge. Remove first face using
766  // edge (and not yet changed)
767 
768  forAll(eFaces, i)
769  {
770  label patchFacei = eFaces[i];
771 
772  label ownerCell = outsideOwner[patchFacei];
773 
774  if (operator[](ownerCell) == meshType)
775  {
776  operator[](ownerCell) = fillType;
777 
778  nChanged++;
779 
780  break;
781  }
782  }
783  }
784  }
785 
786  nTotChanged += nChanged;
787 
788  Pout<< "fillRegionEdges : changed " << nChanged
789  << " cells using multiply connected edges" << endl;
790 
791  if (nChanged == 0)
792  {
793  break;
794  }
795  }
796 
797  return nTotChanged;
798 }
799 
800 
802 (
803  const label meshType,
804  const label fillType,
805  const label maxIter
806 )
807 {
808  label nTotChanged = 0;
809 
810  for (label iter = 0; iter < maxIter; iter++)
811  {
812  // Get interface between meshType cells and non-meshType cells as a list
813  // of faces and for each face the cell which is the meshType.
814  faceList outsideFaces;
815  labelList outsideOwner;
816 
817  getMeshOutside(meshType, outsideFaces, outsideOwner);
818 
819  // Build primitivePatch out of it and check it for problems.
820  primitiveFacePatch fp(outsideFaces, mesh_.points());
821 
822  labelHashSet nonManifoldPoints;
823 
824  // Check for non-manifold points.
825  fp.checkPointManifold(false, &nonManifoldPoints);
826 
827  const Map<label>& meshPointMap = fp.meshPointMap();
828 
829  label nChanged = 0;
830 
831  forAllConstIter(labelHashSet, nonManifoldPoints, iter)
832  {
833  // Find a face on fp using point and remove it.
834  const label patchPointi = meshPointMap[iter.key()];
835 
836  const labelList& pFaces = fp.pointFaces()[patchPointi];
837 
838  // Remove any face using conflicting point. Does first face which
839  // has not yet been done. Could be more intelligent and decide which
840  // one would be best to remove.
841  forAll(pFaces, i)
842  {
843  const label patchFacei = pFaces[i];
844  const label ownerCell = outsideOwner[patchFacei];
845 
846  if (operator[](ownerCell) == meshType)
847  {
848  operator[](ownerCell) = fillType;
849 
850  nChanged++;
851  break;
852  }
853  }
854  }
855 
856  nTotChanged += nChanged;
857 
858  Pout<< "fillRegionPoints : changed " << nChanged
859  << " cells using multiply connected points" << endl;
860 
861  if (nChanged == 0)
862  {
863  break;
864  }
865  }
866 
867  return nTotChanged;
868 }
869 
870 
872 {
873  os << "Cells:" << size() << endl
874  << " notset : " << count(*this, NOTSET) << endl
875  << " cut : " << count(*this, CUT) << endl
876  << " inside : " << count(*this, INSIDE) << endl
877  << " outside : " << count(*this, OUTSIDE) << endl;
878 }
879 
880 
881 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
882 
884 {
886 }
887 
888 
889 // ************************************************************************* //
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
label trimCutCells(const label nLayers, const label meshType, const label fillType)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
static const Form max
Definition: VectorSpace.H:115
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
T & operator[](const label)
Return element of UList.
Definition: UListI.H:167
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
List< face > faceList
Definition: faceListFwd.H:43
label nCells() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Holds information (coordinate and normal) regarding nearest wall point.
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void operator=(const cellClassification &)
fvMesh & mesh
static const Form min
Definition: VectorSpace.H:116
Various functions to operate on Lists.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
cellClassification(const polyMesh &mesh, const meshSearch &meshQuery, const triSurfaceSearch &surfQuery, const pointField &outsidePoints)
Construct from mesh and surface and point(s) on outside.
Helper class to search on triSurface.
A list of faces which address into the list of points.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
&#39;Cuts&#39; a mesh with a surface.
List< edge > edgeList
Definition: edgeList.H:38
label fillRegionPoints(const label meshType, const label fillType, const label maxIter)
Find regionPoints and fill all neighbours. Iterate until nothing.
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const Cmpt & x() const
Definition: VectorI.H:75
const labelListList & pointCells() const
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< small) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &mergedCyclicPolyPatch::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
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
void operator=(const UList< label > &)
Assignment to UList operator. Takes linear time.
Definition: List.C:376
void setSize(const label)
Reset size of List.
Definition: List.C:281
vector point
Point is a vector.
Definition: point.H:41
label growSurface(const label meshType, const label fillType)
Sets vertex neighbours of meshType cells to fillType.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
dimensioned< scalar > mag(const dimensioned< Type > &)
label nPoints() const
label fillRegionEdges(const label meshType, const label fillType, const label maxIter)
Find regionEdges and fill one neighbour. Iterate until nothing.
void writeStats(Ostream &os) const
Write statistics on cell types to Ostream.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
label fillHangingCells(const label meshType, const label fillType, const label maxIter)
Find hanging cells (cells with all points on outside) and set their.
label size() const
Return the number of elements in the UList.
Definition: ListI.H:171
Namespace for OpenFOAM.