boundaryMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "boundaryMesh.H"
30 #include "Time.H"
31 #include "polyMesh.H"
32 #include "repatchPolyTopoChanger.H"
33 #include "faceList.H"
34 #include "indexedOctree.H"
35 #include "treeDataPrimitivePatch.H"
36 #include "triSurface.H"
37 #include "SortableList.H"
38 #include "OFstream.H"
39 #include "primitivePatch.H"
40 #include "indirectPrimitivePatch.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46  defineTypeNameAndDebug(boundaryMesh, 0);
47 }
48 
49 // Normal along which to divide faces into categories (used in getNearest)
50 const Foam::vector Foam::boundaryMesh::splitNormal_(3, 2, 1);
51 
52 // Distance to face tolerance for getNearest
53 const Foam::scalar Foam::boundaryMesh::distanceTol_ = 1e-2;
54 
55 
56 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57 
58 // Returns number of feature edges connected to pointi
59 Foam::label Foam::boundaryMesh::nFeatureEdges(label pointi) const
60 {
61  label nFeats = 0;
62 
63  const labelList& pEdges = mesh().pointEdges()[pointi];
64 
65  forAll(pEdges, pEdgeI)
66  {
67  label edgeI = pEdges[pEdgeI];
68 
69  if (edgeToFeature_[edgeI] != -1)
70  {
71  nFeats++;
72  }
73  }
74  return nFeats;
75 }
76 
77 
78 // Returns next feature edge connected to pointi
79 Foam::label Foam::boundaryMesh::nextFeatureEdge
80 (
81  const label edgeI,
82  const label vertI
83 ) const
84 {
85  const labelList& pEdges = mesh().pointEdges()[vertI];
86 
87  forAll(pEdges, pEdgeI)
88  {
89  label nbrEdgeI = pEdges[pEdgeI];
90 
91  if (nbrEdgeI != edgeI)
92  {
93  label featI = edgeToFeature_[nbrEdgeI];
94 
95  if (featI != -1)
96  {
97  return nbrEdgeI;
98  }
99  }
100  }
101 
102  return -1;
103 }
104 
105 
106 // Finds connected feature edges, starting from startPointi and returns
107 // feature labels (not edge labels). Marks feature edges handled in
108 // featVisited.
109 Foam::labelList Foam::boundaryMesh::collectSegment
110 (
111  const boolList& isFeaturePoint,
112  const label startEdgeI,
113  boolList& featVisited
114 ) const
115 {
116  // Find starting feature point on edge.
117 
118  label edgeI = startEdgeI;
119 
120  const edge& e = mesh().edges()[edgeI];
121 
122  label vertI = e.start();
123 
124  while (!isFeaturePoint[vertI])
125  {
126  // Step to next feature edge
127 
128  edgeI = nextFeatureEdge(edgeI, vertI);
129 
130  if ((edgeI == -1) || (edgeI == startEdgeI))
131  {
132  break;
133  }
134 
135  // Step to next vertex on edge
136 
137  const edge& e = mesh().edges()[edgeI];
138 
139  vertI = e.otherVertex(vertI);
140  }
141 
142  //
143  // Now we have:
144  // edgeI : first edge on this segment
145  // vertI : one of the endpoints of this segment
146  //
147  // Start walking other way and storing edges as we go along.
148  //
149 
150  // Untrimmed storage for current segment
151  labelList featLabels(featureEdges_.size());
152 
153  label featLabelI = 0;
154 
155  label initEdgeI = edgeI;
156 
157  do
158  {
159  // Mark edge as visited
160  label featI = edgeToFeature_[edgeI];
161 
162  if (featI == -1)
163  {
165  << "Problem" << abort(FatalError);
166  }
167  featLabels[featLabelI++] = featI;
168 
169  featVisited[featI] = true;
170 
171  // Step to next vertex on edge
172 
173  const edge& e = mesh().edges()[edgeI];
174 
175  vertI = e.otherVertex(vertI);
176 
177  // Step to next feature edge
178 
179  edgeI = nextFeatureEdge(edgeI, vertI);
180 
181  if ((edgeI == -1) || (edgeI == initEdgeI))
182  {
183  break;
184  }
185  }
186  while (!isFeaturePoint[vertI]);
187 
188 
189  // Trim to size
190  featLabels.setSize(featLabelI);
191 
192  return featLabels;
193 }
194 
195 
196 void Foam::boundaryMesh::markEdges
197 (
198  const label maxDistance,
199  const label edgeI,
200  const label distance,
201  labelList& minDistance,
202  DynamicList<label>& visited
203 ) const
204 {
205  if (distance < maxDistance)
206  {
207  // Don't do anything if reached beyond maxDistance.
208 
209  if (minDistance[edgeI] == -1)
210  {
211  // First visit of edge. Store edge label.
212  visited.append(edgeI);
213  }
214  else if (minDistance[edgeI] <= distance)
215  {
216  // Already done this edge
217  return;
218  }
219 
220  minDistance[edgeI] = distance;
221 
222  const edge& e = mesh().edges()[edgeI];
223 
224  // Do edges connected to e.start
225  const labelList& startEdges = mesh().pointEdges()[e.start()];
226 
227  forAll(startEdges, pEdgeI)
228  {
229  markEdges
230  (
231  maxDistance,
232  startEdges[pEdgeI],
233  distance+1,
234  minDistance,
235  visited
236  );
237  }
238 
239  // Do edges connected to e.end
240  const labelList& endEdges = mesh().pointEdges()[e.end()];
241 
242  forAll(endEdges, pEdgeI)
243  {
244  markEdges
245  (
246  maxDistance,
247  endEdges[pEdgeI],
248  distance+1,
249  minDistance,
250  visited
251  );
252  }
253  }
254 }
255 
256 
257 Foam::label Foam::boundaryMesh::findPatchID
258 (
259  const polyPatchList& patches,
260  const word& patchName
261 ) const
262 {
263  forAll(patches, patchi)
264  {
265  if (patches[patchi].name() == patchName)
266  {
267  return patchi;
268  }
269  }
270 
271  return -1;
272 }
273 
274 
276 {
277  wordList names(patches_.size());
278 
279  forAll(patches_, patchi)
280  {
281  names[patchi] = patches_[patchi].name();
282  }
283  return names;
284 }
285 
286 
287 Foam::label Foam::boundaryMesh::whichPatch
288 (
289  const polyPatchList& patches,
290  const label facei
291 ) const
292 {
293  forAll(patches, patchi)
294  {
295  const polyPatch& pp = patches[patchi];
296 
297  if ((facei >= pp.start()) && (facei < (pp.start() + pp.size())))
298  {
299  return patchi;
300  }
301  }
302  return -1;
303 }
304 
305 
306 // Gets labels of changed faces and propagates them to the edges. Returns
307 // labels of edges changed.
308 Foam::labelList Foam::boundaryMesh::faceToEdge
309 (
310  const boolList& regionEdge,
311  const label region,
312  const labelList& changedFaces,
313  labelList& edgeRegion
314 ) const
315 {
316  labelList changedEdges(mesh().nEdges(), -1);
317  label changedI = 0;
318 
319  forAll(changedFaces, i)
320  {
321  label facei = changedFaces[i];
322 
323  const labelList& fEdges = mesh().faceEdges()[facei];
324 
325  forAll(fEdges, fEdgeI)
326  {
327  label edgeI = fEdges[fEdgeI];
328 
329  if (!regionEdge[edgeI] && (edgeRegion[edgeI] == -1))
330  {
331  edgeRegion[edgeI] = region;
332 
333  changedEdges[changedI++] = edgeI;
334  }
335  }
336  }
337 
338  changedEdges.setSize(changedI);
339 
340  return changedEdges;
341 }
342 
343 
344 // Reverse of faceToEdge: gets edges and returns faces
345 Foam::labelList Foam::boundaryMesh::edgeToFace
346 (
347  const label region,
348  const labelList& changedEdges,
349  labelList& faceRegion
350 ) const
351 {
352  labelList changedFaces(mesh().size(), -1);
353  label changedI = 0;
354 
355  forAll(changedEdges, i)
356  {
357  label edgeI = changedEdges[i];
358 
359  const labelList& eFaces = mesh().edgeFaces()[edgeI];
360 
361  forAll(eFaces, eFacei)
362  {
363  label facei = eFaces[eFacei];
364 
365  if (faceRegion[facei] == -1)
366  {
367  faceRegion[facei] = region;
368 
369  changedFaces[changedI++] = facei;
370  }
371  }
372  }
373 
374  changedFaces.setSize(changedI);
375 
376  return changedFaces;
377 }
378 
379 
380 // Finds area, starting at facei, delimited by borderEdge
381 void Foam::boundaryMesh::markZone
382 (
383  const boolList& borderEdge,
384  label facei,
385  label currentZone,
386  labelList& faceZone
387 ) const
388 {
389  faceZone[facei] = currentZone;
390 
391  // List of faces whose faceZone has been set.
392  labelList changedFaces(1, facei);
393  // List of edges whose faceZone has been set.
394  labelList changedEdges;
395 
396  // Zones on all edges.
397  labelList edgeZone(mesh().nEdges(), -1);
398 
399  while (true)
400  {
401  changedEdges = faceToEdge
402  (
403  borderEdge,
404  currentZone,
405  changedFaces,
406  edgeZone
407  );
408 
409  if (debug)
410  {
411  Pout<< "From changedFaces:" << changedFaces.size()
412  << " to changedEdges:" << changedEdges.size()
413  << endl;
414  }
415 
416  if (changedEdges.empty())
417  {
418  break;
419  }
420 
421  changedFaces = edgeToFace(currentZone, changedEdges, faceZone);
422 
423  if (debug)
424  {
425  Pout<< "From changedEdges:" << changedEdges.size()
426  << " to changedFaces:" << changedFaces.size()
427  << endl;
428  }
429 
430  if (changedFaces.empty())
431  {
432  break;
433  }
434  }
435 }
436 
437 
438 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
439 
441 :
442  meshPtr_(nullptr),
443  patches_(),
444  meshFace_(),
445  featurePoints_(),
446  featureEdges_(),
447  featureToEdge_(),
448  edgeToFeature_(),
449  featureSegments_(),
450  extraEdges_()
451 {}
452 
453 
454 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
457 {
458  meshPtr_.reset(nullptr);
459 }
460 
461 
463 {
464  patches_.clear();
465 
466  patches_.setSize(mesh.boundaryMesh().size());
467 
468  // Number of boundary faces
469  const label nBFaces = mesh.nBoundaryFaces();
470 
471  faceList bFaces(nBFaces);
472 
473  meshFace_.setSize(nBFaces);
474 
475  label bFacei = 0;
476 
477  // Collect all boundary faces.
478  forAll(mesh.boundaryMesh(), patchi)
479  {
480  const polyPatch& pp = mesh.boundaryMesh()[patchi];
481 
482  patches_.set
483  (
484  patchi,
485  new boundaryPatch
486  (
487  pp.name(),
488  patchi,
489  pp.size(),
490  bFacei,
491  pp.type()
492  )
493  );
494 
495  // Collect all faces in global numbering.
496  forAll(pp, patchFacei)
497  {
498  meshFace_[bFacei] = pp.start() + patchFacei;
499 
500  bFaces[bFacei] = pp[patchFacei];
501 
502  bFacei++;
503  }
504  }
505 
506 
507  if (debug)
508  {
509  Pout<< "read : patches now:" << endl;
510 
511  forAll(patches_, patchi)
512  {
513  const boundaryPatch& bp = patches_[patchi];
514 
515  Pout<< " name : " << bp.name() << endl
516  << " size : " << bp.size() << endl
517  << " start : " << bp.start() << endl
518  << " type : " << bp.physicalType() << endl
519  << endl;
520  }
521  }
522 
523  //
524  // Construct single patch for all of boundary
525  //
526 
527  // Temporary primitivePatch to calculate compact points & faces.
528  primitiveFacePatch globalPatch(bFaces, mesh.points());
529 
530  // Store in local(compact) addressing
531  clearOut();
532 
533  meshPtr_.reset
534  (
535  new bMesh(globalPatch.localFaces(), globalPatch.localPoints())
536  );
537 
538  if (debug & 2)
539  {
540  const bMesh& msh = *meshPtr_;
541 
542  Pout<< "** Start of Faces **" << endl;
543 
544  forAll(msh, facei)
545  {
546  const face& f = msh[facei];
547 
548  point ctr(Zero);
549 
550  forAll(f, fp)
551  {
552  ctr += msh.points()[f[fp]];
553  }
554  ctr /= f.size();
555 
556  Pout<< " " << facei
557  << " ctr:" << ctr
558  << " verts:" << f
559  << endl;
560  }
561 
562  Pout<< "** End of Faces **" << endl;
563 
564  Pout<< "** Start of Points **" << endl;
565 
566  forAll(msh.points(), pointi)
567  {
568  Pout<< " " << pointi
569  << " coord:" << msh.points()[pointi]
570  << endl;
571  }
572 
573  Pout<< "** End of Points **" << endl;
574  }
575 
576  // Clear edge storage
577  featurePoints_.clear();
578  featureEdges_.clear();
579 
580  featureToEdge_.clear();
581  edgeToFeature_.resize(meshPtr_->nEdges());
582  edgeToFeature_ = -1;
583 
584  featureSegments_.clear();
585  extraEdges_.clear();
586 }
587 
588 
590 {
591  triSurface surf(fName);
592 
593  if (surf.empty())
594  {
595  return;
596  }
597 
598  // Sort according to region
599  SortableList<label> regions(surf.size());
600 
601  forAll(surf, triI)
602  {
603  regions[triI] = surf[triI].region();
604  }
605  regions.sort();
606 
607  // Determine region mapping.
608  Map<label> regionToBoundaryPatch;
609 
610  label oldRegion = -1111;
611  label boundPatch = 0;
612 
613  forAll(regions, i)
614  {
615  if (regions[i] != oldRegion)
616  {
617  regionToBoundaryPatch.insert(regions[i], boundPatch);
618 
619  oldRegion = regions[i];
620  boundPatch++;
621  }
622  }
623 
624  const geometricSurfacePatchList& surfPatches = surf.patches();
625 
626  patches_.clear();
627 
628  if (surfPatches.size() == regionToBoundaryPatch.size())
629  {
630  // There are as many surface patches as region numbers in triangles
631  // so use the surface patches
632 
633  patches_.setSize(surfPatches.size());
634 
635  // Take over patches, setting size to 0 for now.
636  forAll(surfPatches, patchi)
637  {
638  const geometricSurfacePatch& surfPatch = surfPatches[patchi];
639 
640  patches_.set
641  (
642  patchi,
643  new boundaryPatch
644  (
645  surfPatch.name(),
646  patchi,
647  0,
648  0,
649  surfPatch.geometricType()
650  )
651  );
652  }
653  }
654  else
655  {
656  // There are not enough surface patches. Make up my own.
657 
658  patches_.setSize(regionToBoundaryPatch.size());
659 
660  forAll(patches_, patchi)
661  {
662  patches_.set
663  (
664  patchi,
665  new boundaryPatch
666  (
668  patchi,
669  0,
670  0,
672  )
673  );
674  }
675  }
676 
677  //
678  // Copy according into bFaces according to regions
679  //
680 
681  const labelList& indices = regions.indices();
682 
683  faceList bFaces(surf.size());
684 
685  meshFace_.setSize(surf.size());
686 
687  label bFacei = 0;
688 
689  // Current region number
690  label surfRegion = regions[0];
691  label foamRegion = regionToBoundaryPatch[surfRegion];
692 
693  Pout<< "Surface region " << surfRegion << " becomes boundary patch "
694  << foamRegion << " with name " << patches_[foamRegion].name() << endl;
695 
696 
697  // Index in bFaces of start of current patch
698  label startFacei = 0;
699 
700  forAll(indices, indexI)
701  {
702  label triI = indices[indexI];
703 
704  const labelledTri& tri = surf.localFaces()[triI];
705 
706  if (tri.region() != surfRegion)
707  {
708  // Change of region. We now know the size of the previous one.
709  boundaryPatch& bp = patches_[foamRegion];
710 
711  bp.size() = bFacei - startFacei;
712  bp.start() = startFacei;
713 
714  surfRegion = tri.region();
715  foamRegion = regionToBoundaryPatch[surfRegion];
716 
717  Pout<< "Surface region " << surfRegion << " becomes boundary patch "
718  << foamRegion << " with name " << patches_[foamRegion].name()
719  << endl;
720 
721  startFacei = bFacei;
722  }
723 
724  meshFace_[bFacei] = triI;
725 
726  bFaces[bFacei++] = face(tri);
727  }
728 
729  // Final region
730  boundaryPatch& bp = patches_[foamRegion];
731 
732  bp.size() = bFacei - startFacei;
733  bp.start() = startFacei;
734 
735  //
736  // Construct single primitivePatch for all of boundary
737  //
738 
739  clearOut();
740 
741  // Store compact.
742  meshPtr_.reset(new bMesh(bFaces, surf.localPoints()));
743 
744  // Clear edge storage
745  featurePoints_.clear();
746  featureEdges_.clear();
747 
748  featureToEdge_.clear();
749  edgeToFeature_.resize(meshPtr_->nEdges());
750  edgeToFeature_ = -1;
751 
752  featureSegments_.clear();
753  extraEdges_.clear();
754 }
755 
756 
757 void Foam::boundaryMesh::writeTriSurface(const fileName& fName) const
758 {
759  geometricSurfacePatchList surfPatches(patches_.size());
760 
761  forAll(patches_, patchi)
762  {
763  const boundaryPatch& bp = patches_[patchi];
764 
765  surfPatches[patchi] =
767  (
768  bp.name(),
769  patchi,
770  bp.physicalType()
771  );
772  }
773 
774  //
775  // Simple triangulation.
776  //
777 
778  // Get number of triangles per face
779  labelList nTris(mesh().size());
780 
781  label totalNTris = getNTris(0, mesh().size(), nTris);
782 
783  // Determine per face the starting triangle.
784  labelList startTri(mesh().size());
785 
786  label triI = 0;
787 
788  forAll(mesh(), facei)
789  {
790  startTri[facei] = triI;
791 
792  triI += nTris[facei];
793  }
794 
795  // Triangulate
796  labelList triVerts(3*totalNTris);
797 
798  triangulate(0, mesh().size(), totalNTris, triVerts);
799 
800  // Convert to labelledTri
801 
802  List<labelledTri> tris(totalNTris);
803 
804  triI = 0;
805 
806  forAll(patches_, patchi)
807  {
808  const boundaryPatch& bp = patches_[patchi];
809 
810  forAll(bp, patchFacei)
811  {
812  label facei = bp.start() + patchFacei;
813 
814  label triVertI = 3*startTri[facei];
815 
816  for (label faceTriI = 0; faceTriI < nTris[facei]; faceTriI++)
817  {
818  label v0 = triVerts[triVertI++];
819  label v1 = triVerts[triVertI++];
820  label v2 = triVerts[triVertI++];
821 
822  tris[triI++] = labelledTri(v0, v1, v2, patchi);
823  }
824  }
825  }
826 
827  triSurface surf(tris, surfPatches, mesh().points());
828 
829  OFstream surfStream(fName);
830 
831  surf.write(surfStream);
832 }
833 
834 
835 // Get index in this (boundaryMesh) of face nearest to each boundary face in
836 // pMesh.
837 // Originally all triangles/faces of boundaryMesh would be bunged into
838 // one big octree. Problem was that faces on top of each other, differing
839 // only in sign of normal, could not be found separately. It would always
840 // find only one. We could detect that it was probably finding the wrong one
841 // (based on normal) but could not 'tell' the octree to retrieve the other
842 // one (since they occupy exactly the same space)
843 // So now faces get put into different octrees depending on normal.
844 // !It still will not be possible to differentiate between two faces on top
845 // of each other having the same normal
847 (
848  const primitiveMesh& pMesh,
849  const vector& searchSpan
850 ) const
851 {
852 
853  // Divide faces into two bins acc. to normal
854  // - left of splitNormal
855  // - right ,,
856  DynamicList<label> leftFaces(mesh().size()/2);
857  DynamicList<label> rightFaces(mesh().size()/2);
858 
859  forAll(mesh(), bFacei)
860  {
861  scalar sign = mesh().faceNormals()[bFacei] & splitNormal_;
862 
863  if (sign > -1e-5)
864  {
865  rightFaces.append(bFacei);
866  }
867  if (sign < 1e-5)
868  {
869  leftFaces.append(bFacei);
870  }
871  }
872 
873  leftFaces.shrink();
874  rightFaces.shrink();
875 
876  if (debug)
877  {
878  Pout<< "getNearest :"
879  << " rightBin:" << rightFaces.size()
880  << " leftBin:" << leftFaces.size()
881  << endl;
882  }
883 
884  uindirectPrimitivePatch leftPatch
885  (
886  UIndirectList<face>(mesh(), leftFaces),
887  mesh().points()
888  );
889  uindirectPrimitivePatch rightPatch
890  (
891  UIndirectList<face>(mesh(), rightFaces),
892  mesh().points()
893  );
894 
895 
896  // Overall bb
897  treeBoundBox overallBb(mesh().localPoints());
898 
899  // Extend domain slightly (also makes it 3D if was 2D)
900  // Note asymmetry to avoid having faces align with octree cubes.
901  scalar tol = 1e-6 * overallBb.avgDim();
902 
903  point& bbMin = overallBb.min();
904  bbMin.x() -= tol;
905  bbMin.y() -= tol;
906  bbMin.z() -= tol;
907 
908  point& bbMax = overallBb.max();
909  bbMax.x() += 2*tol;
910  bbMax.y() += 2*tol;
911  bbMax.z() += 2*tol;
912 
913  const scalar planarTol =
915  perturbTol();
916 
917 
918  // Create the octrees
920  <
922  > leftTree
923  (
925  (
926  false, // cacheBb
927  leftPatch,
928  planarTol
929  ),
930  overallBb,
931  10, // maxLevel
932  10, // leafSize
933  3.0 // duplicity
934  );
936  <
938  > rightTree
939  (
941  (
942  false, // cacheBb
943  rightPatch,
944  planarTol
945  ),
946  overallBb,
947  10, // maxLevel
948  10, // leafSize
949  3.0 // duplicity
950  );
951 
952  if (debug)
953  {
954  Pout<< "getNearest : built trees" << endl;
955  }
956 
957 
958  const vectorField& ns = mesh().faceNormals();
959 
960 
961  //
962  // Search nearest triangle centre for every polyMesh boundary face
963  //
964 
965  labelList nearestBFacei(pMesh.nBoundaryFaces());
966 
967  const scalar searchDimSqr = magSqr(searchSpan);
968 
969  forAll(nearestBFacei, patchFacei)
970  {
971  label meshFacei = pMesh.nInternalFaces() + patchFacei;
972 
973  const point& ctr = pMesh.faceCentres()[meshFacei];
974 
975  if (debug && (patchFacei % 1000) == 0)
976  {
977  Pout<< "getNearest : patchFace:" << patchFacei
978  << " meshFacei:" << meshFacei << " ctr:" << ctr << endl;
979  }
980 
981 
982  // Get normal from area vector
983  vector n = pMesh.faceAreas()[meshFacei];
984  scalar area = mag(n);
985  n /= area;
986 
987  scalar typDim = -GREAT;
988  const face& f = pMesh.faces()[meshFacei];
989 
990  forAll(f, fp)
991  {
992  typDim = max(typDim, mag(pMesh.points()[f[fp]] - ctr));
993  }
994 
995  // Search right tree
996  pointIndexHit rightInfo = rightTree.findNearest(ctr, searchDimSqr);
997 
998  // Search left tree. Note: could start from rightDist bounding box
999  // instead of starting from top.
1000  pointIndexHit leftInfo = leftTree.findNearest(ctr, searchDimSqr);
1001 
1002  if (rightInfo.hit())
1003  {
1004  if (leftInfo.hit())
1005  {
1006  // Found in both trees. Compare normals.
1007  label rightFacei = rightFaces[rightInfo.index()];
1008  label leftFacei = leftFaces[leftInfo.index()];
1009 
1010  label rightDist = rightInfo.point().dist(ctr);
1011  label leftDist = leftInfo.point().dist(ctr);
1012 
1013  scalar rightSign = n & ns[rightFacei];
1014  scalar leftSign = n & ns[leftFacei];
1015 
1016  if
1017  (
1018  (rightSign > 0 && leftSign > 0)
1019  || (rightSign < 0 && leftSign < 0)
1020  )
1021  {
1022  // Both same sign. Choose nearest.
1023  if (rightDist < leftDist)
1024  {
1025  nearestBFacei[patchFacei] = rightFacei;
1026  }
1027  else
1028  {
1029  nearestBFacei[patchFacei] = leftFacei;
1030  }
1031  }
1032  else
1033  {
1034  // Differing sign.
1035  // - if both near enough choose one with correct sign
1036  // - otherwise choose nearest.
1037 
1038  // Get local dimension as max of distance between ctr and
1039  // any face vertex.
1040 
1041  typDim *= distanceTol_;
1042 
1043  if (rightDist < typDim && leftDist < typDim)
1044  {
1045  // Different sign and nearby. Choosing matching normal
1046  if (rightSign > 0)
1047  {
1048  nearestBFacei[patchFacei] = rightFacei;
1049  }
1050  else
1051  {
1052  nearestBFacei[patchFacei] = leftFacei;
1053  }
1054  }
1055  else
1056  {
1057  // Different sign but faraway. Choosing nearest.
1058  if (rightDist < leftDist)
1059  {
1060  nearestBFacei[patchFacei] = rightFacei;
1061  }
1062  else
1063  {
1064  nearestBFacei[patchFacei] = leftFacei;
1065  }
1066  }
1067  }
1068  }
1069  else
1070  {
1071  // Found in right but not in left. Choose right regardless if
1072  // correct sign. Note: do we want this?
1073  label rightFacei = rightFaces[rightInfo.index()];
1074  nearestBFacei[patchFacei] = rightFacei;
1075  }
1076  }
1077  else
1078  {
1079  // No face found in right tree.
1080 
1081  if (leftInfo.hit())
1082  {
1083  // Found in left but not in right. Choose left regardless if
1084  // correct sign. Note: do we want this?
1085  nearestBFacei[patchFacei] = leftFaces[leftInfo.index()];
1086  }
1087  else
1088  {
1089  // No face found in left tree.
1090  nearestBFacei[patchFacei] = -1;
1091  }
1092  }
1093  }
1094 
1095  return nearestBFacei;
1096 }
1097 
1098 
1100 (
1101  const labelList& nearest,
1102  const polyBoundaryMesh& oldPatches,
1103  polyMesh& newMesh
1104 ) const
1105 {
1106 
1107  // 2 cases to be handled:
1108  // A- patches in boundaryMesh patches_
1109  // B- patches not in boundaryMesh patches_ but in polyMesh
1110 
1111  // Create maps from patch name to new patch index.
1112  HashTable<label> nameToIndex(2*patches_.size());
1113 
1114  Map<word> indexToName(2*patches_.size());
1115 
1116 
1117  label nNewPatches = patches_.size();
1118 
1119  forAll(oldPatches, oldPatchi)
1120  {
1121  const polyPatch& patch = oldPatches[oldPatchi];
1122  const label newPatchi = findPatchID(patch.name());
1123 
1124  if (newPatchi != -1)
1125  {
1126  nameToIndex.insert(patch.name(), newPatchi);
1127  indexToName.insert(newPatchi, patch.name());
1128  }
1129  }
1130 
1131  // Include all boundaryPatches not yet in nameToIndex (i.e. not in old
1132  // patches)
1133  forAll(patches_, bPatchi)
1134  {
1135  const boundaryPatch& bp = patches_[bPatchi];
1136 
1137  if (!nameToIndex.found(bp.name()))
1138  {
1139  nameToIndex.insert(bp.name(), bPatchi);
1140  indexToName.insert(bPatchi, bp.name());
1141  }
1142  }
1143 
1144  // Pass1:
1145  // Copy names&type of patches (with zero size) from old mesh as far as
1146  // possible. First patch created gets all boundary faces; all others get
1147  // zero faces (repatched later on). Exception is coupled patches which
1148  // keep their size.
1149 
1150  polyPatchList newPatches(nNewPatches);
1151 
1152  label meshFacei = newMesh.nInternalFaces();
1153 
1154  // First patch gets all non-coupled faces
1155  label facesToBeDone = newMesh.nBoundaryFaces();
1156 
1157  forAll(patches_, bPatchi)
1158  {
1159  const boundaryPatch& bp = patches_[bPatchi];
1160 
1161  const label newPatchi = nameToIndex[bp.name()];
1162 
1163  // Find corresponding patch in polyMesh
1164  const label oldPatchi = findPatchID(oldPatches, bp.name());
1165 
1166  if (oldPatchi == -1)
1167  {
1168  // Newly created patch. Gets all or zero faces.
1169  if (debug)
1170  {
1171  Pout<< "patchify : Creating new polyPatch:" << bp.name()
1172  << " type:" << bp.physicalType() << endl;
1173  }
1174 
1175  newPatches.set
1176  (
1177  newPatchi,
1179  (
1180  bp.physicalType(),
1181  bp.name(),
1182  facesToBeDone,
1183  meshFacei,
1184  newPatchi,
1185  newMesh.boundaryMesh()
1186  )
1187  );
1188 
1189  meshFacei += facesToBeDone;
1190 
1191  // first patch gets all boundary faces; all others get 0.
1192  facesToBeDone = 0;
1193  }
1194  else
1195  {
1196  // Existing patch. Gets all or zero faces.
1197  const polyPatch& oldPatch = oldPatches[oldPatchi];
1198 
1199  if (debug)
1200  {
1201  Pout<< "patchify : Cloning existing polyPatch:"
1202  << oldPatch.name() << endl;
1203  }
1204 
1205  newPatches.set
1206  (
1207  newPatchi,
1208  oldPatch.clone
1209  (
1210  newMesh.boundaryMesh(),
1211  newPatchi,
1212  facesToBeDone,
1213  meshFacei
1214  )
1215  );
1216 
1217  meshFacei += facesToBeDone;
1218 
1219  // first patch gets all boundary faces; all others get 0.
1220  facesToBeDone = 0;
1221  }
1222  }
1223 
1224 
1225  if (debug)
1226  {
1227  Pout<< "Patchify : new polyPatch list:" << endl;
1228 
1229  forAll(newPatches, patchi)
1230  {
1231  const polyPatch& newPatch = newPatches[patchi];
1232 
1233  if (debug)
1234  {
1235  Pout<< "polyPatch:" << newPatch.name() << endl
1236  << " type :" << newPatch.typeName << endl
1237  << " size :" << newPatch.size() << endl
1238  << " start:" << newPatch.start() << endl
1239  << " index:" << patchi << endl;
1240  }
1241  }
1242  }
1243 
1244  // Actually add new list of patches
1245  repatchPolyTopoChanger polyMeshRepatcher(newMesh);
1246  polyMeshRepatcher.changePatches(newPatches);
1247 
1248 
1249  // Pass2:
1250  // Change patch type for face
1251 
1252  if (newPatches.size())
1253  {
1254  List<DynamicList<label>> patchFaces(nNewPatches);
1255 
1256  // Give reasonable estimate for size of patches
1257  label nAvgFaces = newMesh.nBoundaryFaces() / nNewPatches;
1258 
1259  forAll(patchFaces, newPatchi)
1260  {
1261  patchFaces[newPatchi].setCapacity(nAvgFaces);
1262  }
1263 
1264  //
1265  // Sort faces acc. to new patch index. Can loop over all old patches
1266  // since will contain all faces.
1267  //
1268 
1269  forAll(oldPatches, oldPatchi)
1270  {
1271  const polyPatch& patch = oldPatches[oldPatchi];
1272 
1273  forAll(patch, patchFacei)
1274  {
1275  // Put face into region given by nearest boundary face
1276 
1277  label meshFacei = patch.start() + patchFacei;
1278 
1279  label bFacei = meshFacei - newMesh.nInternalFaces();
1280 
1281  patchFaces[whichPatch(nearest[bFacei])].append(meshFacei);
1282  }
1283  }
1284 
1285  forAll(patchFaces, newPatchi)
1286  {
1287  patchFaces[newPatchi].shrink();
1288  }
1289 
1290 
1291  // Change patch > 0. (since above we put all faces into the zeroth
1292  // patch)
1293 
1294  for (label newPatchi = 1; newPatchi < patchFaces.size(); newPatchi++)
1295  {
1296  const labelList& pFaces = patchFaces[newPatchi];
1297 
1298  forAll(pFaces, pFacei)
1299  {
1300  polyMeshRepatcher.changePatchID(pFaces[pFacei], newPatchi);
1301  }
1302  }
1303 
1304  polyMeshRepatcher.repatch();
1305  }
1306 }
1307 
1308 
1309 void Foam::boundaryMesh::setFeatureEdges(const scalar minCos)
1310 {
1311  edgeToFeature_.setSize(mesh().nEdges());
1312 
1313  edgeToFeature_ = -1;
1314 
1315  // 1. Mark feature edges
1316 
1317  // Storage for edge labels that are features. Trim later.
1318  featureToEdge_.setSize(mesh().nEdges());
1319 
1320  label featureI = 0;
1321 
1322  if (minCos >= 0.9999)
1323  {
1324  // Select everything
1325  forAll(mesh().edges(), edgeI)
1326  {
1327  edgeToFeature_[edgeI] = featureI;
1328  featureToEdge_[featureI++] = edgeI;
1329  }
1330  }
1331  else
1332  {
1333  forAll(mesh().edges(), edgeI)
1334  {
1335  const labelList& eFaces = mesh().edgeFaces()[edgeI];
1336 
1337  if (eFaces.size() == 2)
1338  {
1339  label face0I = eFaces[0];
1340 
1341  label face1I = eFaces[1];
1342 
1345  //if (whichPatch(face0I) != whichPatch(face1I))
1346  //{
1347  // edgeToFeature_[edgeI] = featureI;
1348  // featureToEdge_[featureI++] = edgeI;
1349  //}
1350  //else
1351  {
1352  const vector& n0 = mesh().faceNormals()[face0I];
1353 
1354  const vector& n1 = mesh().faceNormals()[face1I];
1355 
1356  float cosAng = n0 & n1;
1357 
1358  if (cosAng < minCos)
1359  {
1360  edgeToFeature_[edgeI] = featureI;
1361  featureToEdge_[featureI++] = edgeI;
1362  }
1363  }
1364  }
1365  else
1366  {
1367  //Should not occur: 0 or more than two faces
1368  edgeToFeature_[edgeI] = featureI;
1369  featureToEdge_[featureI++] = edgeI;
1370  }
1371  }
1372  }
1373 
1374  // Trim featureToEdge_ to actual number of edges.
1375  featureToEdge_.setSize(featureI);
1376 
1377  //
1378  // Compact edges i.e. relabel vertices.
1379  //
1380 
1381  featureEdges_.setSize(featureI);
1382  featurePoints_.setSize(mesh().nPoints());
1383 
1384  labelList featToMeshPoint(mesh().nPoints(), -1);
1385 
1386  label featPtI = 0;
1387 
1388  forAll(featureToEdge_, fEdgeI)
1389  {
1390  label edgeI = featureToEdge_[fEdgeI];
1391 
1392  const edge& e = mesh().edges()[edgeI];
1393 
1394  label start = featToMeshPoint[e.start()];
1395 
1396  if (start == -1)
1397  {
1398  featToMeshPoint[e.start()] = featPtI;
1399 
1400  featurePoints_[featPtI] = mesh().points()[e.start()];
1401 
1402  start = featPtI;
1403 
1404  featPtI++;
1405  }
1406 
1407  label end = featToMeshPoint[e.end()];
1408 
1409  if (end == -1)
1410  {
1411  featToMeshPoint[e.end()] = featPtI;
1412 
1413  featurePoints_[featPtI] = mesh().points()[e.end()];
1414 
1415  end = featPtI;
1416 
1417  featPtI++;
1418  }
1419 
1420  // Store with renumbered vertices.
1421  featureEdges_[fEdgeI] = edge(start, end);
1422  }
1423 
1424  // Compact points
1425  featurePoints_.setSize(featPtI);
1426 
1427 
1428  //
1429  // 2. Mark endpoints of feature segments. These are points with
1430  // != 2 feature edges connected.
1431  // Note: can add geometric constraint here as well that if 2 feature
1432  // edges the angle between them should be less than xxx.
1433  //
1434 
1435  boolList isFeaturePoint(mesh().nPoints(), false);
1436 
1437  forAll(featureToEdge_, featI)
1438  {
1439  label edgeI = featureToEdge_[featI];
1440 
1441  const edge& e = mesh().edges()[edgeI];
1442 
1443  if (nFeatureEdges(e.start()) != 2)
1444  {
1445  isFeaturePoint[e.start()] = true;
1446  }
1447 
1448  if (nFeatureEdges(e.end()) != 2)
1449  {
1450  isFeaturePoint[e.end()] = true;
1451  }
1452  }
1453 
1454 
1455  //
1456  // 3: Split feature edges into segments:
1457  // find point with not 2 feature edges -> start of feature segment
1458  //
1459 
1460  DynamicList<labelList> segments;
1461 
1462 
1463  boolList featVisited(featureToEdge_.size(), false);
1464 
1465  do
1466  {
1467  label startFeatI = -1;
1468 
1469  forAll(featVisited, featI)
1470  {
1471  if (!featVisited[featI])
1472  {
1473  startFeatI = featI;
1474 
1475  break;
1476  }
1477  }
1478 
1479  if (startFeatI == -1)
1480  {
1481  // No feature lines left.
1482  break;
1483  }
1484 
1485  segments.append
1486  (
1487  collectSegment
1488  (
1489  isFeaturePoint,
1490  featureToEdge_[startFeatI],
1491  featVisited
1492  )
1493  );
1494  }
1495  while (true);
1496 
1497 
1498  //
1499  // Store in *this
1500  //
1501  featureSegments_.setSize(segments.size());
1502 
1503  forAll(featureSegments_, segmentI)
1504  {
1505  featureSegments_[segmentI] = segments[segmentI];
1506  }
1507 }
1508 
1509 
1510 void Foam::boundaryMesh::setExtraEdges(const label edgeI)
1511 {
1512  labelList minDistance(mesh().nEdges(), -1);
1513 
1514  // All edge labels encountered
1515  DynamicList<label> visitedEdges;
1516 
1517  // Floodfill from edgeI starting from distance 0. Stop at distance.
1518  markEdges(8, edgeI, 0, minDistance, visitedEdges);
1519 
1520  // Set edge labels to display
1521  extraEdges_.transfer(visitedEdges);
1522 }
1523 
1524 
1525 Foam::label Foam::boundaryMesh::whichPatch(const label facei) const
1526 {
1527  forAll(patches_, patchi)
1528  {
1529  const boundaryPatch& bp = patches_[patchi];
1530 
1531  if ((facei >= bp.start()) && (facei < (bp.start() + bp.size())))
1532  {
1533  return patchi;
1534  }
1535  }
1536 
1538  << "Cannot find face " << facei << " in list of boundaryPatches "
1539  << patches_
1540  << abort(FatalError);
1541 
1542  return -1;
1543 }
1544 
1545 
1546 Foam::label Foam::boundaryMesh::findPatchID(const word& patchName) const
1547 {
1548  forAll(patches_, patchi)
1549  {
1550  if (patches_[patchi].name() == patchName)
1551  {
1552  return patchi;
1553  }
1554  }
1555 
1556  return -1;
1557 }
1558 
1559 
1560 void Foam::boundaryMesh::addPatch(const word& patchName)
1561 {
1562  patches_.setSize(patches_.size() + 1);
1563 
1564  // Add empty patch at end of patch list.
1565 
1566  label patchi = patches_.size()-1;
1567 
1568  boundaryPatch* bpPtr = new boundaryPatch
1569  (
1570  patchName,
1571  patchi,
1572  0,
1573  mesh().size(),
1574  "empty"
1575  );
1576 
1577  patches_.set(patchi, bpPtr);
1578 
1579  if (debug)
1580  {
1581  Pout<< "addPatch : patches now:" << endl;
1582 
1583  forAll(patches_, patchi)
1584  {
1585  const boundaryPatch& bp = patches_[patchi];
1586 
1587  Pout<< " name : " << bp.name() << endl
1588  << " size : " << bp.size() << endl
1589  << " start : " << bp.start() << endl
1590  << " type : " << bp.physicalType() << endl
1591  << endl;
1592  }
1593  }
1594 }
1595 
1596 
1597 void Foam::boundaryMesh::deletePatch(const word& patchName)
1598 {
1599  const label delPatchi = findPatchID(patchName);
1600 
1601  if (delPatchi == -1)
1602  {
1604  << "Can't find patch named " << patchName
1605  << abort(FatalError);
1606  }
1607 
1608  if (patches_[delPatchi].size())
1609  {
1611  << "Trying to delete non-empty patch " << patchName
1612  << endl << "Current size:" << patches_[delPatchi].size()
1613  << abort(FatalError);
1614  }
1615 
1616  PtrList<boundaryPatch> newPatches(patches_.size() - 1);
1617 
1618  for (label patchi = 0; patchi < delPatchi; patchi++)
1619  {
1620  newPatches.set(patchi, patches_[patchi].clone());
1621  }
1622 
1623  // Move patches down, starting from delPatchi.
1624 
1625  for (label patchi = delPatchi + 1; patchi < patches_.size(); patchi++)
1626  {
1627  newPatches.set(patchi - 1, patches_[patchi].clone());
1628  }
1629 
1630  patches_.clear();
1631 
1632  patches_ = newPatches;
1633 
1634  if (debug)
1635  {
1636  Pout<< "deletePatch : patches now:" << endl;
1637 
1638  forAll(patches_, patchi)
1639  {
1640  const boundaryPatch& bp = patches_[patchi];
1641 
1642  Pout<< " name : " << bp.name() << endl
1643  << " size : " << bp.size() << endl
1644  << " start : " << bp.start() << endl
1645  << " type : " << bp.physicalType() << endl
1646  << endl;
1647  }
1648  }
1649 }
1650 
1651 
1653 (
1654  const word& patchName,
1655  const word& patchType
1656 )
1657 {
1658  const label changeI = findPatchID(patchName);
1659 
1660  if (changeI == -1)
1661  {
1663  << "Can't find patch named " << patchName
1664  << abort(FatalError);
1665  }
1666 
1667 
1668  // Cause we can't reassign to individual PtrList elems ;-(
1669  // work on copy
1670 
1671 
1672  PtrList<boundaryPatch> newPatches(patches_.size());
1673 
1674  forAll(patches_, patchi)
1675  {
1676  if (patchi == changeI)
1677  {
1678  // Create copy but for type
1679  const boundaryPatch& oldBp = patches_[patchi];
1680 
1681  boundaryPatch* bpPtr = new boundaryPatch
1682  (
1683  oldBp.name(),
1684  oldBp.index(),
1685  oldBp.size(),
1686  oldBp.start(),
1687  patchType
1688  );
1689 
1690  newPatches.set(patchi, bpPtr);
1691  }
1692  else
1693  {
1694  // Create copy
1695  newPatches.set(patchi, patches_[patchi].clone());
1696  }
1697  }
1698 
1699  patches_ = newPatches;
1700 }
1701 
1702 
1704 (
1705  const labelList& patchIDs,
1706  labelList& oldToNew
1707 )
1708 {
1709  if (patchIDs.size() != mesh().size())
1710  {
1712  << "List of patchIDs not equal to number of faces." << endl
1713  << "PatchIDs size:" << patchIDs.size()
1714  << " nFaces::" << mesh().size()
1715  << abort(FatalError);
1716  }
1717 
1718  // Count number of faces for each patch
1719 
1720  labelList nFaces(patches_.size(), Zero);
1721 
1722  forAll(patchIDs, facei)
1723  {
1724  label patchID = patchIDs[facei];
1725 
1726  if (patchID < 0 || patchID >= patches_.size())
1727  {
1729  << "PatchID " << patchID << " out of range"
1730  << abort(FatalError);
1731  }
1732  nFaces[patchID]++;
1733  }
1734 
1735 
1736  // Determine position in faces_ for each patch
1737 
1738  labelList startFace(patches_.size());
1739 
1740  startFace[0] = 0;
1741 
1742  for (label patchi = 1; patchi < patches_.size(); patchi++)
1743  {
1744  startFace[patchi] = startFace[patchi-1] + nFaces[patchi-1];
1745  }
1746 
1747  // Update patch info
1748  PtrList<boundaryPatch> newPatches(patches_.size());
1749 
1750  forAll(patches_, patchi)
1751  {
1752  const boundaryPatch& bp = patches_[patchi];
1753 
1754  newPatches.set
1755  (
1756  patchi,
1757  new boundaryPatch
1758  (
1759  bp.name(),
1760  patchi,
1761  nFaces[patchi],
1762  startFace[patchi],
1763  bp.physicalType()
1764  )
1765  );
1766  }
1767  patches_ = newPatches;
1768 
1769  if (debug)
1770  {
1771  Pout<< "changeFaces : patches now:" << endl;
1772 
1773  forAll(patches_, patchi)
1774  {
1775  const boundaryPatch& bp = patches_[patchi];
1776 
1777  Pout<< " name : " << bp.name() << endl
1778  << " size : " << bp.size() << endl
1779  << " start : " << bp.start() << endl
1780  << " type : " << bp.physicalType() << endl
1781  << endl;
1782  }
1783  }
1784 
1785 
1786  // Construct face mapping array
1787  oldToNew.setSize(patchIDs.size());
1788 
1789  forAll(patchIDs, facei)
1790  {
1791  int patchID = patchIDs[facei];
1792 
1793  oldToNew[facei] = startFace[patchID]++;
1794  }
1795 
1796  // Copy faces into correct position and maintain label of original face
1797  faceList newFaces(mesh().size());
1798 
1799  labelList newMeshFace(mesh().size());
1800 
1801  forAll(oldToNew, facei)
1802  {
1803  newFaces[oldToNew[facei]] = mesh()[facei];
1804  newMeshFace[oldToNew[facei]] = meshFace_[facei];
1805  }
1806 
1807  // Reconstruct 'mesh' from new faces and (copy of) existing points.
1808 
1809  std::unique_ptr<bMesh> newMeshPtr(new bMesh(newFaces, mesh().points()));
1810 
1811  // Reset meshFace_ to new ordering.
1812  meshFace_.transfer(newMeshFace);
1813 
1814 
1815  // Remove old PrimitivePatch on meshPtr_.
1816  clearOut();
1817 
1818  // And insert new 'mesh'.
1819  meshPtr_ = std::move(newMeshPtr);
1820 }
1821 
1822 
1823 Foam::label Foam::boundaryMesh::getNTris(const label facei) const
1824 {
1825  const face& f = mesh()[facei];
1826 
1827  return f.nTriangles(mesh().points());
1828 }
1829 
1830 
1831 Foam::label Foam::boundaryMesh::getNTris
1832 (
1833  const label startFacei,
1834  const label nFaces,
1835  labelList& nTris
1836 ) const
1837 {
1838  label totalNTris = 0;
1839 
1840  nTris.setSize(nFaces);
1841 
1842  for (label i = 0; i < nFaces; i++)
1843  {
1844  label faceNTris = getNTris(startFacei + i);
1845 
1846  nTris[i] = faceNTris;
1847 
1848  totalNTris += faceNTris;
1849  }
1850  return totalNTris;
1851 }
1852 
1853 
1854 // Simple triangulation of face subset. Stores vertices in tris[] as three
1855 // consecutive vertices per triangle.
1857 (
1858  const label startFacei,
1859  const label nFaces,
1860  const label totalNTris,
1861  labelList& triVerts
1862 ) const
1863 {
1864  // Triangulate faces.
1865  triVerts.setSize(3*totalNTris);
1866 
1867  label vertI = 0;
1868 
1869  for (label i = 0; i < nFaces; i++)
1870  {
1871  label facei = startFacei + i;
1872 
1873  const face& f = mesh()[facei];
1874 
1875  // Have face triangulate itself (results in faceList)
1876  faceList triFaces(f.nTriangles(mesh().points()));
1877 
1878  label nTri = 0;
1879 
1880  f.triangles(mesh().points(), nTri, triFaces);
1881 
1882  // Copy into triVerts
1883 
1884  forAll(triFaces, triFacei)
1885  {
1886  const face& triF = triFaces[triFacei];
1887 
1888  triVerts[vertI++] = triF[0];
1889  triVerts[vertI++] = triF[1];
1890  triVerts[vertI++] = triF[2];
1891  }
1892  }
1893 }
1894 
1895 
1896 // Number of local points in subset.
1898 (
1899  const label startFacei,
1900  const label nFaces
1901 ) const
1902 {
1904  (
1905  SubList<face>(mesh(), nFaces, startFacei),
1906  mesh().points()
1907  );
1909  return patch.nPoints();
1910 }
1911 
1912 
1913 // Triangulation of face subset in local coords.
1915 (
1916  const label startFacei,
1917  const label nFaces,
1918  const label totalNTris,
1919  labelList& triVerts,
1920  labelList& localToGlobal
1921 ) const
1922 {
1924  (
1925  SubList<face>(mesh(), nFaces, startFacei),
1926  mesh().points()
1927  );
1928 
1929  // Map from local to mesh().points() addressing
1930  localToGlobal = patch.meshPoints();
1931 
1932  // Triangulate patch faces.
1933  triVerts.setSize(3*totalNTris);
1934 
1935  label vertI = 0;
1936 
1937  for (label i = 0; i < nFaces; i++)
1938  {
1939  // Local face
1940  const face& f = patch.localFaces()[i];
1941 
1942  // Have face triangulate itself (results in faceList)
1943  faceList triFaces(f.nTriangles(patch.localPoints()));
1944 
1945  label nTri = 0;
1946 
1947  f.triangles(patch.localPoints(), nTri, triFaces);
1948 
1949  // Copy into triVerts
1950 
1951  forAll(triFaces, triFacei)
1952  {
1953  const face& triF = triFaces[triFacei];
1954 
1955  triVerts[vertI++] = triF[0];
1956  triVerts[vertI++] = triF[1];
1957  triVerts[vertI++] = triF[2];
1958  }
1959  }
1960 }
1961 
1962 
1964 (
1965  const labelList& protectedEdges,
1966  const label seedFacei,
1967  boolList& visited
1968 ) const
1969 {
1970  boolList protectedEdge(mesh().nEdges(), false);
1971 
1972  forAll(protectedEdges, i)
1973  {
1974  protectedEdge[protectedEdges[i]] = true;
1975  }
1976 
1977 
1978  // Initialize zone for all faces to -1
1979  labelList currentZone(mesh().size(), -1);
1980 
1981  // Mark with 0 all faces reachable from seedFacei
1982  markZone(protectedEdge, seedFacei, 0, currentZone);
1983 
1984  // Set in visited all reached ones.
1985  visited.setSize(mesh().size());
1986 
1987  forAll(currentZone, facei)
1988  {
1989  if (currentZone[facei] == 0)
1990  {
1991  visited[facei] = true;
1992  }
1993  else
1994  {
1995  visited[facei] = false;
1996  }
1997  }
1998 }
1999 
2000 
2001 // ************************************************************************* //
label size() const
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
dimensionedScalar sign(const dimensionedScalar &ds)
const labelListList & pointEdges() const
Return point-edge addressing.
virtual const fileName & name() const
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
void addPatch(const word &patchName)
Add to back of patch list.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
List of names generated by calling name() for each list item and filtered for matches.
A class for handling file names.
Definition: fileName.H:71
void read(const polyMesh &)
Read from boundaryMesh of polyMesh.
Definition: boundaryMesh.C:455
void setFeatureEdges(const scalar minCos)
Set featureEdges, edgeToFeature, featureSegments according.
Like polyPatch but without reference to mesh. Used in boundaryMesh to hold data on patches...
Definition: boundaryPatch.H:50
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const Field< point_type > & localPoints() const
Return pointField of points in patch.
label getNPoints(const label startFacei, const label nFaces) const
Number of points used in face subset.
List< geometricSurfacePatch > geometricSurfacePatchList
A List of geometricSurfacePatch.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
const labelListList & faceEdges() const
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
const labelListList & pointEdges() const
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
void triangulate(const label startFacei, const label nFaces, const label totalNTris, labelList &triVerts) const
Simple triangulation of face subset. TotalNTris is total number.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:420
List< face > faceList
A List of faces.
Definition: faceListFwd.H:41
label getNTris(const label facei) const
Simple triangulation of face subset. Returns number of triangles.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
void writeTriSurface(const fileName &) const
Write to file.
Definition: boundaryMesh.C:750
virtual const pointField & points() const =0
Return mesh points.
void patchify(const labelList &nearest, const polyBoundaryMesh &oldPatches, polyMesh &newMesh) const
Take over patches onto polyMesh from nearest face in *this.
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:45
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
labelList getNearest(const primitiveMesh &pMesh, const vector &searchSpan) const
Get bMesh index of nearest face for every boundary face in.
Definition: boundaryMesh.C:840
void markFaces(const labelList &protectedEdges, const label facei, boolList &visited) const
const point & min() const noexcept
Minimum describing the bounding box.
Definition: boundBoxI.H:155
label start() const
void write(Ostream &os) const
Write to Ostream in simple OpenFOAM format.
Definition: triSurfaceIO.C:331
void changePatchType(const word &patchName, const word &type)
Change patch.
scalar avgDim() const
Average length/height/width dimension.
Definition: boundBoxI.H:221
void deletePatch(const word &patchName)
Delete from patch list.
boundaryMesh()
Default construct.
Definition: boundaryMesh.C:433
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1066
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
static word defaultName(const label n=-1)
Default patch name: "patch" or "patchN".
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition: polyPatch.H:56
const point & max() const noexcept
Maximum describing the bounding box.
Definition: boundBoxI.H:161
A list of faces which address into the list of points.
const point_type & point() const noexcept
Return point, no checks.
A List obtained as a section of another List.
Definition: SubList.H:50
const geometricSurfacePatchList & patches() const noexcept
Definition: triSurface.H:510
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:289
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
const pointField & points
const wordList area
Standard area field types (scalar, vector, tensor, etc)
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:109
A class for handling words, derived from Foam::string.
Definition: word.H:63
label nPoints
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:99
Encapsulation of data needed to search on PrimitivePatches.
label nInternalFaces() const noexcept
Number of internal faces.
Vector< scalar > vector
Definition: vector.H:57
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:558
label index() const noexcept
Return the hit index.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:163
wordList patchNames() const
Get names of patches.
Definition: boundaryMesh.C:268
const word & name() const noexcept
The patch name.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:427
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:193
const word & physicalType() const noexcept
The (optional) physical type of the patch.
defineTypeNameAndDebug(combustionModel, 0)
static constexpr const char *const emptyType
The name for an &#39;empty&#39; type.
labelList f(nPoints)
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const vectorField & faceCentres() const
void setExtraEdges(const label edgeI)
Set extraEdges to edges &#39;near&#39; to edgeI. Uses point-edge walk.
List< word > wordList
A List of words.
Definition: fileName.H:58
bool hit() const noexcept
Is there a hit?
vector point
Point is a vector.
Definition: point.H:37
Non-pointer based hierarchical recursive searching.
void readTriSurface(const fileName &)
Read from triSurface.
Definition: boundaryMesh.C:582
const vectorField & faceAreas() const
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:57
virtual const faceList & faces() const =0
Return faces.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:441
const polyBoundaryMesh & patches
const std::string patch
OpenFOAM patch number as a std::string.
void triangulateLocal(const label startFacei, const label nFaces, const label totalNTris, labelList &triVerts, labelList &localToGlobal) const
Same as triangulate but in local vertex numbering.
Identifies a surface patch/zone by name and index, with geometric type.
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
label n
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
Definition: polyPatchNew.C:28
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
const bMesh & mesh() const
Definition: boundaryMesh.H:259
List< label > labelList
A List of labels.
Definition: List.H:62
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Triangulated surface description with patch information.
Definition: triSurface.H:72
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points) ...
Definition: bMesh.H:42
List< bool > boolList
A List of bools.
Definition: List.H:60
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const labelListList & edgeFaces() const
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
void changeFaces(const labelList &patchIDs, labelList &oldToNew)
Recalculate face ordering and patches. Return old to new.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157