extendedEdgeMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2021 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 "extendedEdgeMesh.H"
30 #include "surfaceFeatures.H"
31 #include "triSurface.H"
32 #include "Random.H"
33 #include "Time.H"
34 #include "OBJstream.H"
35 #include "DynamicField.H"
36 #include "edgeMeshFormatsCore.H"
37 #include "IOmanip.H"
38 #include "searchableSurface.H"
39 #include "triSurfaceMesh.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45  defineTypeNameAndDebug(extendedEdgeMesh, 0);
46 }
47 
48 
49 const Foam::Enum
50 <
52 >
54 ({
55  { pointStatus::CONVEX, "convex" },
56  { pointStatus::CONCAVE, "concave" },
57  { pointStatus::MIXED, "mixed" },
58  { pointStatus::NONFEATURE, "nonFeature" },
59 });
60 
61 
62 const Foam::Enum
63 <
65 >
67 ({
68  { edgeStatus::EXTERNAL, "external" },
69  { edgeStatus::INTERNAL, "internal" },
70  { edgeStatus::FLAT, "flat" },
71  { edgeStatus::OPEN, "open" },
72  { edgeStatus::MULTIPLE, "multiple" },
73  { edgeStatus::NONE, "none" },
74 });
75 
76 
77 const Foam::Enum
78 <
80 >
82 ({
83  { sideVolumeType::INSIDE, "inside" },
84  { sideVolumeType::OUTSIDE, "outside" },
85  { sideVolumeType::BOTH, "both" },
86  { sideVolumeType::NEITHER, "neither" },
87 });
88 
89 
91  Foam::cos(degToRad(0.1));
92 
93 
95 
97 
98 
99 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
100 
102 {
103  return wordHashSet(*fileExtensionConstructorTablePtr_);
104 }
105 
106 
108 {
109  return wordHashSet(*writefileExtensionMemberFunctionTablePtr_);
110 }
111 
112 
113 bool Foam::extendedEdgeMesh::canReadType(const word& fileType, bool verbose)
114 {
115  return edgeMeshFormatsCore::checkSupport
116  (
117  readTypes(),
118  fileType,
119  verbose,
120  "reading"
121  );
122 }
123 
124 
125 bool Foam::extendedEdgeMesh::canWriteType(const word& fileType, bool verbose)
126 {
127  return edgeMeshFormatsCore::checkSupport
128  (
129  writeTypes(),
130  fileType,
131  verbose,
132  "writing"
133  );
134 }
135 
136 
138 {
139  word ext(name.ext());
140  if (ext == "gz")
141  {
142  ext = name.lessExt().ext();
143  }
144  return canReadType(ext, verbose);
145 }
146 
147 
148 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
149 
152 (
153  label ptI
154 ) const
155 {
156  const labelList& ptEds(pointEdges()[ptI]);
157 
158  label nPtEds = ptEds.size();
159  label nExternal = 0;
160  label nInternal = 0;
161 
162  if (nPtEds == 0)
163  {
164  // There are no edges attached to the point, this is a problem
165  return NONFEATURE;
166  }
167 
168  forAll(ptEds, i)
169  {
170  edgeStatus edStat = getEdgeStatus(ptEds[i]);
171 
172  if (edStat == EXTERNAL)
173  {
174  nExternal++;
175  }
176  else if (edStat == INTERNAL)
177  {
178  nInternal++;
179  }
180  }
181 
182  if (nExternal == nPtEds)
183  {
184  return CONVEX;
185  }
186  else if (nInternal == nPtEds)
187  {
188  return CONCAVE;
189  }
190 
191  return MIXED;
192 }
193 
194 
196 (
197  const searchableSurface& surf,
198 
199  labelList& pointMap,
200  labelList& edgeMap,
201  labelList& pointsFromEdge,
202  labelList& oldEdge,
203  labelList& surfTri
204 )
205 {
206  const edgeList& edges = this->edges();
207  const pointField& points = this->points();
208 
209 
210  List<List<pointIndexHit>> edgeHits(edges.size());
211  {
212  pointField start(edges.size());
213  pointField end(edges.size());
214  forAll(edges, edgeI)
215  {
216  const edge& e = edges[edgeI];
217  start[edgeI] = points[e[0]];
218  end[edgeI] = points[e[1]];
219  }
220  surf.findLineAll(start, end, edgeHits);
221  }
222 
223  // Count number of hits
224  label nHits = 0;
225  forAll(edgeHits, edgeI)
226  {
227  nHits += edgeHits[edgeI].size();
228  }
229 
230  DynamicField<point> newPoints(points);
231  DynamicList<label> newToOldPoint(identity(points.size()));
232 
233  newPoints.setCapacity(newPoints.size()+nHits);
234  newToOldPoint.setCapacity(newPoints.capacity());
235 
236  DynamicList<edge> newEdges(edges);
237  DynamicList<label> newToOldEdge(identity(edges.size()));
238 
239  newEdges.setCapacity(newEdges.size()+nHits);
240  newToOldEdge.setCapacity(newEdges.capacity());
241 
242  // Information on additional points
243  DynamicList<label> dynPointsFromEdge(nHits);
244  DynamicList<label> dynOldEdge(nHits);
245  DynamicList<label> dynSurfTri(nHits);
246 
247  forAll(edgeHits, edgeI)
248  {
249  const List<pointIndexHit>& eHits = edgeHits[edgeI];
250 
251  if (eHits.size())
252  {
253  label prevPtI = edges[edgeI][0];
254  forAll(eHits, eHitI)
255  {
256  label newPtI = newPoints.size();
257 
258  newPoints.append(eHits[eHitI].hitPoint());
259  newToOldPoint.append(edges[edgeI][0]); // map from start point
260  dynPointsFromEdge.append(newPtI);
261  dynOldEdge.append(edgeI);
262  dynSurfTri.append(eHits[eHitI].index());
263 
264  if (eHitI == 0)
265  {
266  newEdges[edgeI] = edge(prevPtI, newPtI);
267  }
268  else
269  {
270  newEdges.append(edge(prevPtI, newPtI));
271  newToOldEdge.append(edgeI);
272  }
273  prevPtI = newPtI;
274  }
275  newEdges.append(edge(prevPtI, edges[edgeI][1]));
276  newToOldEdge.append(edgeI);
277  }
278  }
279 
281  allPoints.transfer(newPoints);
282  pointMap.transfer(newToOldPoint);
283 
284  edgeList allEdges;
285  allEdges.transfer(newEdges);
286  edgeMap.transfer(newToOldEdge);
287 
288  pointsFromEdge.transfer(dynPointsFromEdge);
289  oldEdge.transfer(dynOldEdge);
290  surfTri.transfer(dynSurfTri);
291 
292  // Update local information
293  autoMap(allPoints, allEdges, pointMap, edgeMap);
294 }
295 
296 
298 (
299  const searchableSurface& surf,
300  const volumeType volType, // side to keep
301  labelList& pointMap, // sub to old points
302  labelList& edgeMap // sub to old edges
303 )
304 {
305  const edgeList& edges = this->edges();
306  const pointField& points = this->points();
307 
308  // Test edge centres for inside/outside
309  if (volType == volumeType::INSIDE || volType == volumeType::OUTSIDE)
310  {
311  pointField edgeCentres(edges.size());
312  forAll(edgeCentres, edgeI)
313  {
314  const edge& e = edges[edgeI];
315  edgeCentres[edgeI] = e.centre(points);
316  }
317  List<volumeType> volTypes;
318  surf.getVolumeType(edgeCentres, volTypes);
319 
320  // Extract edges on correct side
321  edgeMap.setSize(edges.size());
322  label compactEdgeI = 0;
323 
324  forAll(volTypes, edgeI)
325  {
326  if (volTypes[edgeI] == volType)
327  {
328  edgeMap[compactEdgeI++] = edgeI;
329  }
330  }
331  edgeMap.setSize(compactEdgeI);
332 
333  // Extract used points
334  labelList pointToCompact(points.size(), -1);
335  forAll(edgeMap, i)
336  {
337  const edge& e = edges[edgeMap[i]];
338  pointToCompact[e[0]] = labelMax; // tag with a value
339  pointToCompact[e[1]] = labelMax;
340  }
341 
342  pointMap.setSize(points.size());
343  label compactPointI = 0;
344  forAll(pointToCompact, pointI)
345  {
346  if (pointToCompact[pointI] != -1)
347  {
348  pointToCompact[pointI] = compactPointI;
349  pointMap[compactPointI++] = pointI;
350  }
351  }
352  pointMap.setSize(compactPointI);
353  pointField subPoints(points, pointMap);
354 
355  // Renumber edges
356  edgeList subEdges(edgeMap.size());
357  forAll(edgeMap, i)
358  {
359  const edge& e = edges[edgeMap[i]];
360  subEdges[i][0] = pointToCompact[e[0]];
361  subEdges[i][1] = pointToCompact[e[1]];
362  }
363 
364  // Reset primitives and map other quantities
365  autoMap(subPoints, subEdges, pointMap, edgeMap);
366  }
367  else
368  {
369  pointMap = identity(points.size());
370  edgeMap = identity(edges.size());
371  }
372 }
373 
374 
375 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
376 
378 :
379  edgeMesh(),
380  concaveStart_(0),
381  mixedStart_(0),
382  nonFeatureStart_(0),
383  internalStart_(0),
384  flatStart_(0),
385  openStart_(0),
386  multipleStart_(0),
387  normals_(0),
388  normalVolumeTypes_(0),
389  edgeDirections_(0),
390  normalDirections_(0),
391  edgeNormals_(0),
392  featurePointNormals_(0),
393  featurePointEdges_(0),
394  regionEdges_(0),
395  pointTree_(nullptr),
396  edgeTree_(nullptr),
397  edgeTreesByType_()
398 {}
399 
400 
402 :
403  edgeMesh(),
404  concaveStart_(-1),
405  mixedStart_(-1),
406  nonFeatureStart_(-1),
407  internalStart_(-1),
408  flatStart_(-1),
409  openStart_(-1),
410  multipleStart_(-1),
411  normals_(0),
412  normalVolumeTypes_(0),
413  edgeDirections_(0),
414  normalDirections_(0),
415  edgeNormals_(0),
416  featurePointNormals_(0),
417  featurePointEdges_(0),
418  regionEdges_(0),
419  pointTree_(nullptr),
420  edgeTree_(nullptr),
421  edgeTreesByType_()
422 {}
423 
424 
426 :
427  edgeMesh(fem),
428  concaveStart_(fem.concaveStart()),
429  mixedStart_(fem.mixedStart()),
430  nonFeatureStart_(fem.nonFeatureStart()),
431  internalStart_(fem.internalStart()),
432  flatStart_(fem.flatStart()),
433  openStart_(fem.openStart()),
434  multipleStart_(fem.multipleStart()),
435  normals_(fem.normals()),
436  normalVolumeTypes_(fem.normalVolumeTypes()),
437  edgeDirections_(fem.edgeDirections()),
438  normalDirections_(fem.normalDirections()),
439  edgeNormals_(fem.edgeNormals()),
440  featurePointNormals_(fem.featurePointNormals()),
441  featurePointEdges_(fem.featurePointEdges()),
442  regionEdges_(fem.regionEdges()),
443  pointTree_(nullptr),
444  edgeTree_(nullptr),
445  edgeTreesByType_()
446 {}
447 
448 
450 {
451  is >> *this;
452 }
453 
454 
456 (
457  const pointField& points,
458  const edgeList& edges
459 )
460 :
462 {
463  this->storedPoints() = points;
464  this->storedEdges() = edges;
465 }
466 
467 
469 (
470  pointField&& points,
471  edgeList&& edges
472 )
473 :
475 {
476  this->storedPoints().transfer(points);
477  this->storedEdges().transfer(edges);
478 }
479 
480 
482 (
483  const surfaceFeatures& sFeat,
484  const boolList& surfBaffleRegions
485 )
486 :
488 {
489  // Extract and reorder the data from surfaceFeatures
490  const triSurface& surf = sFeat.surface();
491  const labelList& featureEdges = sFeat.featureEdges();
492  const labelList& featurePoints = sFeat.featurePoints();
493 
494  // Get a labelList of all the featureEdges that are region edges
495  const labelList regionFeatureEdges(identity(sFeat.nRegionEdges()));
496 
497  sortPointsAndEdges
498  (
499  surf,
500  featureEdges,
501  regionFeatureEdges,
502  featurePoints
503  );
504 
505  const labelListList& edgeFaces = surf.edgeFaces();
506 
507  normalVolumeTypes_.setSize(normals_.size());
508 
509  // Noting when the normal of a face has been used so not to duplicate
510  labelList faceMap(surf.size(), -1);
511 
512  label nAdded = 0;
513 
514  forAll(featureEdges, i)
515  {
516  label sFEI = featureEdges[i];
517 
518  // Pick up the faces adjacent to the feature edge
519  const labelList& eFaces = edgeFaces[sFEI];
520 
521  forAll(eFaces, j)
522  {
523  label eFI = eFaces[j];
524 
525  // Check to see if the points have been already used
526  if (faceMap[eFI] == -1)
527  {
528  normalVolumeTypes_[nAdded++] =
529  (
530  surfBaffleRegions[surf[eFI].region()]
531  ? BOTH
532  : INSIDE
533  );
534 
535  faceMap[eFI] = nAdded - 1;
536  }
537  }
538  }
539 }
540 
541 
543 (
545  const labelUList& featureEdges,
546  const labelUList& regionFeatureEdges,
547  const labelUList& featurePoints
548 )
549 :
551 {
552  sortPointsAndEdges
553  (
554  surf,
555  featureEdges,
556  regionFeatureEdges,
557  featurePoints
558  );
559 }
560 
561 
563 (
564  const pointField& pts,
565  const edgeList& eds,
566  label concaveStart,
567  label mixedStart,
568  label nonFeatureStart,
569  label internalStart,
570  label flatStart,
571  label openStart,
572  label multipleStart,
573  const vectorField& normals,
574  const List<sideVolumeType>& normalVolumeTypes,
575  const vectorField& edgeDirections,
576  const labelListList& normalDirections,
577  const labelListList& edgeNormals,
578  const labelListList& featurePointNormals,
579  const labelListList& featurePointEdges,
580  const labelList& regionEdges
581 )
582 :
583  edgeMesh(pts, eds),
584  concaveStart_(concaveStart),
585  mixedStart_(mixedStart),
586  nonFeatureStart_(nonFeatureStart),
587  internalStart_(internalStart),
588  flatStart_(flatStart),
589  openStart_(openStart),
590  multipleStart_(multipleStart),
591  normals_(normals),
592  normalVolumeTypes_(normalVolumeTypes),
593  edgeDirections_(edgeDirections),
594  normalDirections_(normalDirections),
595  edgeNormals_(edgeNormals),
596  featurePointNormals_(featurePointNormals),
597  featurePointEdges_(featurePointEdges),
598  regionEdges_(regionEdges),
599  pointTree_(nullptr),
600  edgeTree_(nullptr),
601  edgeTreesByType_()
602 {}
603 
604 
606 (
607  const fileName& name,
608  const word& fileType
609 )
610 :
612 {
613  read(name, fileType);
614 }
615 
616 
618 :
620 {
621  read(name);
622 }
623 
624 
625 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
626 
628 {
629  word ext(name.ext());
630  if (ext == "gz")
631  {
632  fileName unzipName = name.lessExt();
633  return read(unzipName, unzipName.ext());
634  }
635 
636  return read(name, ext);
637 }
638 
639 
641 (
642  const fileName& name,
643  const word& fileType
644 )
645 {
646  // Read via selector mechanism
647  transfer(*New(name, fileType));
648  return true;
649 }
650 
651 
653 (
654  const point& sample,
655  scalar searchDistSqr,
656  pointIndexHit& info
657 ) const
658 {
659  info = pointTree().findNearest
660  (
661  sample,
662  searchDistSqr
663  );
664 }
665 
666 
668 (
669  const point& sample,
670  scalar searchDistSqr,
671  pointIndexHit& info
672 ) const
673 {
674  info = edgeTree().findNearest
675  (
676  sample,
677  searchDistSqr
678  );
679 }
680 
681 
683 (
684  const pointField& samples,
685  const scalarField& searchDistSqr,
686  List<pointIndexHit>& info
687 ) const
688 {
689  info.setSize(samples.size());
690 
691  forAll(samples, i)
692  {
693  nearestFeatureEdge
694  (
695  samples[i],
696  searchDistSqr[i],
697  info[i]
698  );
699  }
700 }
701 
702 
704 (
705  const point& sample,
706  const scalarField& searchDistSqr,
707  List<pointIndexHit>& info
708 ) const
709 {
710  const PtrList<indexedOctree<treeDataEdge>>& edgeTrees = edgeTreesByType();
711 
712  info.setSize(edgeTrees.size());
713 
714  labelList sliceStarts(edgeTrees.size());
715 
716  sliceStarts[0] = externalStart_;
717  sliceStarts[1] = internalStart_;
718  sliceStarts[2] = flatStart_;
719  sliceStarts[3] = openStart_;
720  sliceStarts[4] = multipleStart_;
721 
722  forAll(edgeTrees, i)
723  {
724  info[i] = edgeTrees[i].findNearest
725  (
726  sample,
727  searchDistSqr[i]
728  );
729 
730  // The index returned by the indexedOctree is local to the slice of
731  // edges it was supplied with, return the index to the value in the
732  // complete edge list
733 
734  info[i].setIndex(info[i].index() + sliceStarts[i]);
735  }
736 }
737 
738 
740 (
741  const point& sample,
742  scalar searchRadiusSqr,
743  List<pointIndexHit>& info
744 ) const
745 {
746  // Pick up all the feature points that intersect the search sphere
747  labelList elems = pointTree().findSphere
748  (
749  sample,
750  searchRadiusSqr
751  );
752 
753  DynamicList<pointIndexHit> dynPointHit(elems.size());
754 
755  forAll(elems, elemI)
756  {
757  label index = elems[elemI];
758  label ptI = pointTree().shapes().pointLabels()[index];
759  const point& pt = points()[ptI];
760 
761  pointIndexHit nearHit(true, pt, index);
762 
763  dynPointHit.append(nearHit);
764  }
765 
766  info.transfer(dynPointHit);
767 }
768 
769 
771 (
772  const point& sample,
773  const scalar searchRadiusSqr,
774  List<pointIndexHit>& info
775 ) const
776 {
777  const PtrList<indexedOctree<treeDataEdge>>& edgeTrees = edgeTreesByType();
778 
779  info.setSize(edgeTrees.size());
780 
781  labelList sliceStarts(edgeTrees.size());
782 
783  sliceStarts[0] = externalStart_;
784  sliceStarts[1] = internalStart_;
785  sliceStarts[2] = flatStart_;
786  sliceStarts[3] = openStart_;
787  sliceStarts[4] = multipleStart_;
788 
789  DynamicList<pointIndexHit> dynEdgeHit(edgeTrees.size()*3);
790 
791  // Loop over all the feature edge types
792  forAll(edgeTrees, i)
793  {
794  // Pick up all the edges that intersect the search sphere
795  labelList elems = edgeTrees[i].findSphere
796  (
797  sample,
798  searchRadiusSqr
799  );
800 
801  forAll(elems, elemI)
802  {
803  label index = elems[elemI];
804  label edgeI = edgeTrees[i].shapes().edgeLabels()[index];
805  const edge& e = edges()[edgeI];
806 
807  pointHit hitPoint = e.line(points()).nearestDist(sample);
808 
809  label hitIndex = index + sliceStarts[i];
810 
811  pointIndexHit nearHit
812  (
813  hitPoint.hit(),
814  hitPoint.rawPoint(),
815  hitIndex
816  );
817 
818  dynEdgeHit.append(nearHit);
819  }
820  }
821 
822  info.transfer(dynEdgeHit);
823 }
824 
825 
828 {
829  if (!pointTree_)
830  {
831  Random rndGen(17301893);
832 
833  // Slightly extended bb. Slightly off-centred just so on symmetric
834  // geometry there are less face/edge aligned items.
835  treeBoundBox bb
836  (
837  treeBoundBox(points()).extend(rndGen, 1e-4)
838  );
839 
840  bb.min() -= point::uniform(ROOTVSMALL);
841  bb.max() += point::uniform(ROOTVSMALL);
842 
843  const labelList featurePointLabels = identity(nonFeatureStart_);
844 
845  pointTree_.reset
846  (
848  (
850  (
851  points(),
852  featurePointLabels
853  ),
854  bb, // bb
855  8, // maxLevel
856  10, // leafsize
857  3.0 // duplicity
858  )
859  );
860  }
861 
862  return *pointTree_;
863 }
864 
865 
868 {
869  if (!edgeTree_)
870  {
871  Random rndGen(17301893);
872 
873  // Slightly extended bb. Slightly off-centred just so on symmetric
874  // geometry there are less face/edge aligned items.
875  treeBoundBox bb
876  (
877  treeBoundBox(points()).extend(rndGen, 1e-4)
878  );
879 
880  bb.min() -= point::uniform(ROOTVSMALL);
881  bb.max() += point::uniform(ROOTVSMALL);
882 
883  labelList allEdges(identity(edges().size()));
884 
885  edgeTree_.reset
886  (
888  (
890  (
891  false, // cachebb
892  edges(), // edges
893  points(), // points
894  allEdges // selected edges
895  ),
896  bb, // bb
897  8, // maxLevel
898  10, // leafsize
899  3.0 // duplicity
900  )
901  );
902  }
903 
904  return *edgeTree_;
905 }
906 
907 
910 {
911  if (edgeTreesByType_.empty())
912  {
913  Random rndGen(872141);
914 
915  // Slightly extended bb. Slightly off-centred just so on symmetric
916  // geometry there are less face/edge aligned items.
917  treeBoundBox bb
918  (
919  treeBoundBox(points()).extend(rndGen, 1e-4)
920  );
921 
922  bb.min() -= point::uniform(ROOTVSMALL);
923  bb.max() += point::uniform(ROOTVSMALL);
924 
925  labelListList sliceEdges(nEdgeTypes);
926 
927  // External edges
928  sliceEdges[0] =
929  identity((internalStart_ - externalStart_), externalStart_);
930 
931  // Internal edges
932  sliceEdges[1] = identity((flatStart_ - internalStart_), internalStart_);
933 
934  // Flat edges
935  sliceEdges[2] = identity((openStart_ - flatStart_), flatStart_);
936 
937  // Open edges
938  sliceEdges[3] = identity((multipleStart_ - openStart_), openStart_);
939 
940  // Multiple edges
941  sliceEdges[4] =
942  identity((edges().size() - multipleStart_), multipleStart_);
943 
944 
945  edgeTreesByType_.resize(nEdgeTypes);
946 
947  forAll(edgeTreesByType_, i)
948  {
949  edgeTreesByType_.set
950  (
951  i,
953  (
955  (
956  false, // cachebb
957  edges(), // edges
958  points(), // points
959  sliceEdges[i] // selected edges
960  ),
961  bb, // bb
962  8, // maxLevel
963  10, // leafsize
964  3.0 // duplicity
965  )
966  );
967  }
968  }
969 
970  return edgeTreesByType_;
971 }
972 
973 
975 {
976  if (&mesh == this)
977  {
978  return; // Self-transfer is a no-op
979  }
980 
982 
983  concaveStart_ = mesh.concaveStart_;
984  mixedStart_ = mesh.mixedStart_;
985  nonFeatureStart_ = mesh.nonFeatureStart_;
986  internalStart_ = mesh.internalStart_;
987  flatStart_ = mesh.flatStart_;
988  openStart_ = mesh.openStart_;
989  multipleStart_ = mesh.multipleStart_;
990  normals_.transfer(mesh.normals_);
991  normalVolumeTypes_.transfer(mesh.normalVolumeTypes_);
992  edgeDirections_.transfer(mesh.edgeDirections_);
993  normalDirections_.transfer(mesh.normalDirections_);
994  edgeNormals_.transfer(mesh.edgeNormals_);
995  featurePointNormals_.transfer(mesh.featurePointNormals_);
996  featurePointEdges_.transfer(mesh.featurePointEdges_);
997  regionEdges_.transfer(mesh.regionEdges_);
998  pointTree_ = std::move(mesh.pointTree_);
999  edgeTree_ = std::move(mesh.edgeTree_);
1000  edgeTreesByType_.transfer(mesh.edgeTreesByType_);
1001 
1002  mesh.clear();
1003 }
1004 
1005 
1007 {
1008  edgeMesh::clear();
1009  concaveStart_ = 0;
1010  mixedStart_ = 0;
1011  nonFeatureStart_ = 0;
1012  internalStart_ = 0;
1013  flatStart_ = 0;
1014  openStart_ = 0;
1015  multipleStart_ = 0;
1016  normals_.clear();
1017  normalVolumeTypes_.clear();
1018  edgeDirections_.clear();
1019  normalDirections_.clear();
1020  edgeNormals_.clear();
1021  featurePointNormals_.clear();
1022  featurePointEdges_.clear();
1023  regionEdges_.clear();
1024  pointTree_.reset(nullptr);
1025  edgeTree_.reset(nullptr);
1026  edgeTreesByType_.clear();
1027 }
1028 
1029 
1031 {
1032  // Points
1033  // ~~~~~~
1034 
1035  // From current points into combined points
1036  labelList reversePointMap(points().size());
1037  labelList reverseFemPointMap(fem.points().size());
1038 
1039  label newPointi = 0;
1040  for (label i = 0; i < concaveStart(); i++)
1041  {
1042  reversePointMap[i] = newPointi++;
1043  }
1044  for (label i = 0; i < fem.concaveStart(); i++)
1045  {
1046  reverseFemPointMap[i] = newPointi++;
1047  }
1048 
1049  // Concave
1050  label newConcaveStart = newPointi;
1051  for (label i = concaveStart(); i < mixedStart(); i++)
1052  {
1053  reversePointMap[i] = newPointi++;
1054  }
1055  for (label i = fem.concaveStart(); i < fem.mixedStart(); i++)
1056  {
1057  reverseFemPointMap[i] = newPointi++;
1058  }
1059 
1060  // Mixed
1061  label newMixedStart = newPointi;
1062  for (label i = mixedStart(); i < nonFeatureStart(); i++)
1063  {
1064  reversePointMap[i] = newPointi++;
1065  }
1066  for (label i = fem.mixedStart(); i < fem.nonFeatureStart(); i++)
1067  {
1068  reverseFemPointMap[i] = newPointi++;
1069  }
1070 
1071  // Non-feature
1072  label newNonFeatureStart = newPointi;
1073  for (label i = nonFeatureStart(); i < points().size(); i++)
1074  {
1075  reversePointMap[i] = newPointi++;
1076  }
1077  for (label i = fem.nonFeatureStart(); i < fem.points().size(); i++)
1078  {
1079  reverseFemPointMap[i] = newPointi++;
1080  }
1081 
1082  pointField newPoints(newPointi);
1083  newPoints.rmap(points(), reversePointMap);
1084  newPoints.rmap(fem.points(), reverseFemPointMap);
1085 
1086 
1087  // Edges
1088  // ~~~~~
1089 
1090  // From current edges into combined edges
1091  labelList reverseEdgeMap(edges().size());
1092  labelList reverseFemEdgeMap(fem.edges().size());
1093 
1094  // External
1095  label newEdgeI = 0;
1096  for (label i = 0; i < internalStart(); i++)
1097  {
1098  reverseEdgeMap[i] = newEdgeI++;
1099  }
1100  for (label i = 0; i < fem.internalStart(); i++)
1101  {
1102  reverseFemEdgeMap[i] = newEdgeI++;
1103  }
1104 
1105  // Internal
1106  label newInternalStart = newEdgeI;
1107  for (label i = internalStart(); i < flatStart(); i++)
1108  {
1109  reverseEdgeMap[i] = newEdgeI++;
1110  }
1111  for (label i = fem.internalStart(); i < fem.flatStart(); i++)
1112  {
1113  reverseFemEdgeMap[i] = newEdgeI++;
1114  }
1115 
1116  // Flat
1117  label newFlatStart = newEdgeI;
1118  for (label i = flatStart(); i < openStart(); i++)
1119  {
1120  reverseEdgeMap[i] = newEdgeI++;
1121  }
1122  for (label i = fem.flatStart(); i < fem.openStart(); i++)
1123  {
1124  reverseFemEdgeMap[i] = newEdgeI++;
1125  }
1126 
1127  // Open
1128  label newOpenStart = newEdgeI;
1129  for (label i = openStart(); i < multipleStart(); i++)
1130  {
1131  reverseEdgeMap[i] = newEdgeI++;
1132  }
1133  for (label i = fem.openStart(); i < fem.multipleStart(); i++)
1134  {
1135  reverseFemEdgeMap[i] = newEdgeI++;
1136  }
1137 
1138  // Multiple
1139  label newMultipleStart = newEdgeI;
1140  for (label i = multipleStart(); i < edges().size(); i++)
1141  {
1142  reverseEdgeMap[i] = newEdgeI++;
1143  }
1144  for (label i = fem.multipleStart(); i < fem.edges().size(); i++)
1145  {
1146  reverseFemEdgeMap[i] = newEdgeI++;
1147  }
1148 
1149  edgeList newEdges(newEdgeI);
1150  forAll(edges(), i)
1151  {
1152  const edge& e = edges()[i];
1153  newEdges[reverseEdgeMap[i]] = edge
1154  (
1155  reversePointMap[e[0]],
1156  reversePointMap[e[1]]
1157  );
1158  }
1159  forAll(fem.edges(), i)
1160  {
1161  const edge& e = fem.edges()[i];
1162  newEdges[reverseFemEdgeMap[i]] = edge
1163  (
1164  reverseFemPointMap[e[0]],
1165  reverseFemPointMap[e[1]]
1166  );
1167  }
1168 
1169  pointField newEdgeDirections
1170  (
1171  edgeDirections().size()
1172  + fem.edgeDirections().size()
1173  );
1174  newEdgeDirections.rmap(edgeDirections(), reverseEdgeMap);
1175  newEdgeDirections.rmap(fem.edgeDirections(), reverseFemEdgeMap);
1176 
1177 
1178  // Normals
1179  // ~~~~~~~
1180 
1181  // Combine normals
1182  DynamicField<point> newNormals
1183  (
1184  normals().size()
1185  + fem.normals().size()
1186  );
1187  newNormals.append(normals());
1188  newNormals.append(fem.normals());
1189 
1190 
1191  // Combine and re-index into newNormals
1192  labelListList newEdgeNormals
1193  (
1194  edgeNormals().size()
1195  + fem.edgeNormals().size()
1196  );
1197 
1199  (
1200  newEdgeNormals,
1201  SubList<label>(reverseEdgeMap, edgeNormals().size())
1202  ) = edgeNormals();
1204  (
1205  newEdgeNormals,
1206  SubList<label>(reverseFemEdgeMap, fem.edgeNormals().size())
1207  ) = fem.edgeNormals();
1208 
1209  forAll(fem.edgeNormals(), i)
1210  {
1211  const label mapI = reverseFemEdgeMap[i];
1212  labelList& en = newEdgeNormals[mapI];
1213  forAll(en, j)
1214  {
1215  en[j] += normals().size();
1216  }
1217  }
1218 
1219 
1220  // Combine and re-index into newFeaturePointNormals
1221  labelListList newFeaturePointNormals
1222  (
1223  featurePointNormals().size()
1224  + fem.featurePointNormals().size()
1225  );
1226 
1227  // Note: featurePointNormals only go up to nonFeatureStart
1229  (
1230  newFeaturePointNormals,
1231  SubList<label>(reversePointMap, featurePointNormals().size())
1232  ) = featurePointNormals();
1234  (
1235  newFeaturePointNormals,
1236  SubList<label>(reverseFemPointMap, fem.featurePointNormals().size())
1237  ) = fem.featurePointNormals();
1238 
1239  forAll(fem.featurePointNormals(), i)
1240  {
1241  const label mapI = reverseFemPointMap[i];
1242  labelList& fn = newFeaturePointNormals[mapI];
1243  forAll(fn, j)
1244  {
1245  fn[j] += normals().size();
1246  }
1247  }
1248 
1249 
1250  // Combine regionEdges
1251  DynamicList<label> newRegionEdges
1252  (
1253  regionEdges().size()
1254  + fem.regionEdges().size()
1255  );
1256  forAll(regionEdges(), i)
1257  {
1258  newRegionEdges.append(reverseEdgeMap[regionEdges()[i]]);
1259  }
1260  forAll(fem.regionEdges(), i)
1261  {
1262  newRegionEdges.append(reverseFemEdgeMap[fem.regionEdges()[i]]);
1263  }
1264 
1265 
1266  // Assign
1267  // ~~~~~~
1268 
1269  // Transfer
1270  concaveStart_ = newConcaveStart;
1271  mixedStart_ = newMixedStart;
1272  nonFeatureStart_ = newNonFeatureStart;
1273 
1274  // Reset points and edges
1275  {
1276  edgeMesh newmesh(std::move(newPoints), std::move(newEdges));
1277  edgeMesh::transfer(newmesh);
1278  }
1279 
1280  // Transfer
1281  internalStart_ = newInternalStart;
1282  flatStart_ = newFlatStart;
1283  openStart_ = newOpenStart;
1284  multipleStart_ = newMultipleStart;
1285 
1286  edgeDirections_.transfer(newEdgeDirections);
1287 
1288  normals_.transfer(newNormals);
1289  edgeNormals_.transfer(newEdgeNormals);
1290  featurePointNormals_.transfer(newFeaturePointNormals);
1291 
1292  regionEdges_.transfer(newRegionEdges);
1293 
1294  pointTree_.reset(nullptr);
1295  edgeTree_.reset(nullptr);
1296  edgeTreesByType_.clear();
1297 }
1298 
1299 
1301 {
1302  // Points
1303  // ~~~~~~
1304 
1305  // From current points into new points
1306  labelList reversePointMap(identity(points().size()));
1307 
1308  // Flip convex and concave points
1309 
1310  label newPointi = 0;
1311  // Concave points become convex
1312  for (label i = concaveStart(); i < mixedStart(); i++)
1313  {
1314  reversePointMap[i] = newPointi++;
1315  }
1316  // Convex points become concave
1317  label newConcaveStart = newPointi;
1318  for (label i = 0; i < concaveStart(); i++)
1319  {
1320  reversePointMap[i] = newPointi++;
1321  }
1322 
1323 
1324  // Edges
1325  // ~~~~~~
1326 
1327  // From current edges into new edges
1328  labelList reverseEdgeMap(identity(edges().size()));
1329 
1330  // Flip external and internal edges
1331 
1332  label newEdgeI = 0;
1333  // Internal become external
1334  for (label i = internalStart(); i < flatStart(); i++)
1335  {
1336  reverseEdgeMap[i] = newEdgeI++;
1337  }
1338  // External become internal
1339  label newInternalStart = newEdgeI;
1340  for (label i = 0; i < internalStart(); i++)
1341  {
1342  reverseEdgeMap[i] = newEdgeI++;
1343  }
1344 
1345 
1346  pointField newPoints(points().size());
1347  newPoints.rmap(points(), reversePointMap);
1348 
1349  edgeList newEdges(edges().size());
1350  forAll(edges(), i)
1351  {
1352  const edge& e = edges()[i];
1353  newEdges[reverseEdgeMap[i]] = edge
1354  (
1355  reversePointMap[e[0]],
1356  reversePointMap[e[1]]
1357  );
1358  }
1359 
1360 
1361  // Normals are flipped
1362  // ~~~~~~~~~~~~~~~~~~~
1363 
1364  pointField newEdgeDirections(edges().size());
1365  newEdgeDirections.rmap(-1.0*edgeDirections(), reverseEdgeMap);
1366 
1367  pointField newNormals(-1.0*normals());
1368 
1369  labelListList newEdgeNormals(edgeNormals().size());
1370  UIndirectList<labelList>(newEdgeNormals, reverseEdgeMap) = edgeNormals();
1371 
1372  labelListList newFeaturePointNormals(featurePointNormals().size());
1373 
1374  // Note: featurePointNormals only go up to nonFeatureStart
1376  (
1377  newFeaturePointNormals,
1378  SubList<label>(reversePointMap, featurePointNormals().size())
1379  ) = featurePointNormals();
1380 
1381  labelList newRegionEdges(regionEdges().size());
1382  forAll(regionEdges(), i)
1383  {
1384  newRegionEdges[i] = reverseEdgeMap[regionEdges()[i]];
1385  }
1386 
1387  // Transfer
1388  concaveStart_ = newConcaveStart;
1389 
1390  // Reset points and edges
1391  {
1392  edgeMesh newmesh(std::move(newPoints), std::move(newEdges));
1393  edgeMesh::transfer(newmesh);
1394  }
1395 
1396  // Transfer
1397  internalStart_ = newInternalStart;
1398 
1399  edgeDirections_.transfer(newEdgeDirections);
1400  normals_.transfer(newNormals);
1401  edgeNormals_.transfer(newEdgeNormals);
1402  featurePointNormals_.transfer(newFeaturePointNormals);
1403  regionEdges_.transfer(newRegionEdges);
1404 
1405  pointTree_.reset(nullptr);
1406  edgeTree_.reset(nullptr);
1407  edgeTreesByType_.clear();
1408 }
1409 
1410 
1413  const pointField& subPoints,
1414  const edgeList& subEdges,
1415  const labelList& pointMap,
1416  const labelList& edgeMap
1417 )
1418 {
1419  // Determine slicing for subEdges
1420  label subIntStart = edgeMap.size();
1421  label subFlatStart = edgeMap.size();
1422  label subOpenStart = edgeMap.size();
1423  label subMultipleStart = edgeMap.size();
1424 
1425  forAll(edgeMap, subEdgeI)
1426  {
1427  label edgeI = edgeMap[subEdgeI];
1428  if (edgeI >= internalStart() && subIntStart == edgeMap.size())
1429  {
1430  subIntStart = subEdgeI;
1431  }
1432  if (edgeI >= flatStart() && subFlatStart == edgeMap.size())
1433  {
1434  subFlatStart = subEdgeI;
1435  }
1436  if (edgeI >= openStart() && subOpenStart == edgeMap.size())
1437  {
1438  subOpenStart = subEdgeI;
1439  }
1440  if (edgeI >= multipleStart() && subMultipleStart == edgeMap.size())
1441  {
1442  subMultipleStart = subEdgeI;
1443  }
1444  }
1445 
1446 
1447  // Determine slicing for subPoints
1448 
1449  label subConcaveStart = pointMap.size();
1450  label subMixedStart = pointMap.size();
1451  label subNonFeatStart = pointMap.size();
1452 
1453  forAll(pointMap, subPointI)
1454  {
1455  label pointI = pointMap[subPointI];
1456  if (pointI >= concaveStart() && subConcaveStart == pointMap.size())
1457  {
1458  subConcaveStart = subPointI;
1459  }
1460  if (pointI >= mixedStart() && subMixedStart == pointMap.size())
1461  {
1462  subMixedStart = subPointI;
1463  }
1464  if
1465  (
1466  pointI >= nonFeatureStart()
1467  && subNonFeatStart == pointMap.size()
1468  )
1469  {
1470  subNonFeatStart = subPointI;
1471  }
1472  }
1473 
1474 
1475 
1476  // Compact region edges
1477  labelList subRegionEdges;
1478  {
1479  bitSet isRegionEdge(edges().size(), regionEdges());
1480 
1481  DynamicList<label> newRegionEdges(regionEdges().size());
1482  forAll(edgeMap, subEdgeI)
1483  {
1484  if (isRegionEdge.test(edgeMap[subEdgeI]))
1485  {
1486  newRegionEdges.append(subEdgeI);
1487  }
1488  }
1489  subRegionEdges.transfer(newRegionEdges);
1490  }
1491 
1492 
1493  labelListList subFeaturePointEdges;
1494  if (featurePointEdges().size())
1495  {
1496  subFeaturePointEdges.setSize(subNonFeatStart);
1497  for (label subPointI = 0; subPointI < subNonFeatStart; subPointI++)
1498  {
1499  label pointI = pointMap[subPointI];
1500  const labelList& pEdges = featurePointEdges()[pointI];
1501 
1502  labelList& subPEdges = subFeaturePointEdges[subPointI];
1503  subPEdges.setSize(pEdges.size());
1504 
1505  if (pEdges.size())
1506  {
1507  forAll(pEdges, i)
1508  {
1509  subPEdges[i] = edgeMap[pEdges[i]];
1510  }
1511  }
1512  }
1513  }
1514 
1515 
1516  vectorField subEdgeDirections(edgeDirections(), edgeMap);
1517 
1518  // Find used normals
1519  labelList reverseNormalMap(normals().size(), -1);
1520  DynamicList<label> normalMap(normals().size());
1521 
1522  {
1523  bitSet isSubNormal(normals().size());
1524  for (label subPointI = 0; subPointI < subNonFeatStart; subPointI++)
1525  {
1526  label pointI = pointMap[subPointI];
1527  const labelList& pNormals = featurePointNormals()[pointI];
1528 
1529  isSubNormal.set(pNormals);
1530  }
1531  forAll(edgeMap, subEdgeI)
1532  {
1533  label edgeI = edgeMap[subEdgeI];
1534  const labelList& eNormals = edgeNormals()[edgeI];
1535 
1536  isSubNormal.set(eNormals);
1537  }
1538 
1539  forAll(isSubNormal, normalI)
1540  {
1541  if (isSubNormal.test(normalI))
1542  {
1543  label subNormalI = normalMap.size();
1544  reverseNormalMap[normalI] = subNormalI;
1545  normalMap.append(subNormalI);
1546  }
1547  }
1548  }
1549 
1550 
1551  // Use compaction map on data referencing normals
1552  labelListList subNormalDirections;
1553 
1554  if (normalDirections().size())
1555  {
1556  subNormalDirections.setSize(edgeMap.size());
1557 
1558  forAll(edgeMap, subEdgeI)
1559  {
1560  label edgeI = edgeMap[subEdgeI];
1561  const labelList& eNormals = normalDirections()[edgeI];
1562 
1563  labelList& subNormals = subNormalDirections[subEdgeI];
1564  subNormals.setSize(eNormals.size());
1565  forAll(eNormals, i)
1566  {
1567  if (eNormals[i] >= 0)
1568  {
1569  subNormals[i] = reverseNormalMap[eNormals[i]];
1570  }
1571  else
1572  {
1573  subNormals[i] = -1;
1574  }
1575  }
1576  }
1577  }
1578 
1579  labelListList subEdgeNormals(edgeMap.size());
1580  forAll(edgeMap, subEdgeI)
1581  {
1582  label edgeI = edgeMap[subEdgeI];
1583  const labelList& eNormals = edgeNormals()[edgeI];
1584  labelList& subNormals = subEdgeNormals[subEdgeI];
1585 
1586  subNormals = labelUIndList(reverseNormalMap, eNormals);
1587  }
1588 
1589  labelListList subPointNormals(pointMap.size());
1590  for (label subPointI = 0; subPointI < subNonFeatStart; subPointI++)
1591  {
1592  label pointI = pointMap[subPointI];
1593  const labelList& pNormals = featurePointNormals()[pointI];
1594  labelList& subNormals = subPointNormals[subPointI];
1595 
1596  subNormals = labelUIndList(reverseNormalMap, pNormals);
1597  }
1598 
1599  // Use compaction map to compact normal data
1600  vectorField subNormals(normals(), normalMap);
1601 
1602  List<extendedEdgeMesh::sideVolumeType> subNormalVolumeTypes;
1603  if (normalVolumeTypes().size())
1604  {
1605  subNormalVolumeTypes =
1607  (
1608  normalVolumeTypes(),
1609  normalMap
1610  );
1611  }
1612 
1613  extendedEdgeMesh subMesh
1614  (
1615  subPoints,
1616  subEdges,
1617 
1618  // Feature points slices
1619  subConcaveStart,
1620  subMixedStart,
1621  subNonFeatStart,
1622 
1623  // Feature edges slices
1624  subIntStart,
1625  subFlatStart,
1626  subOpenStart,
1627  subMultipleStart,
1628 
1629  // All normals
1630  subNormals,
1631  subNormalVolumeTypes,
1632 
1633  // Per edge edge vector
1634  subEdgeDirections,
1635 
1636  // Per edge list of normal indices
1637  subNormalDirections,
1638  // Per edge list of normal indices
1639  subEdgeNormals,
1640 
1641  // Per point list of normal indices
1642  subPointNormals,
1643  subFeaturePointEdges,
1644  subRegionEdges
1645  );
1646 
1647  transfer(subMesh);
1648 }
1649 
1650 
1653  const searchableSurface& surf,
1654  const volumeType volType,
1655  labelList& pointMap,
1656  labelList& edgeMap
1657 )
1658 {
1659  // Cut edges with the other surfaces
1660 
1661  labelList allPointMap; // from all to original point
1662  labelList allEdgeMap; // from all to original edge
1663 
1664  labelList pointsFromEdge; // list of new points created by cutting
1665  labelList oldEdge; // for each of these points the original edge
1666  labelList surfTri; // for each of these points the surface triangle
1667  cut
1668  (
1669  surf,
1670  allPointMap,
1671  allEdgeMap,
1672  pointsFromEdge,
1673  oldEdge,
1674  surfTri
1675  );
1676 
1677  const label nOldPoints = points().size();
1678 
1679  // Remove outside edges and compact
1680 
1681  labelList subPointMap; // sub to old points
1682  labelList subEdgeMap; // sub to old edges
1683  select(surf, volType, subPointMap, subEdgeMap);
1684 
1685  // Update overall point maps
1686  pointMap = labelUIndList(allPointMap, subPointMap);
1687  edgeMap = labelUIndList(allEdgeMap, subEdgeMap);
1688 
1689  // Extract current point and edge status
1690  List<edgeStatus> edgeStat(edges().size());
1691  List<pointStatus> pointStat(points().size());
1692  forAll(edgeStat, edgeI)
1693  {
1694  edgeStat[edgeI] = getEdgeStatus(edgeI);
1695  }
1696  forAll(pointStat, pointI)
1697  {
1698  pointStat[pointI] = getPointStatus(pointI);
1699  }
1700 
1701  // Re-classify exposed points (from cutting)
1702  labelList oldPointToIndex(nOldPoints, -1);
1703  forAll(pointsFromEdge, i)
1704  {
1705  oldPointToIndex[pointsFromEdge[i]] = i;
1706  }
1707  forAll(subPointMap, pointI)
1708  {
1709  label oldPointI = subPointMap[pointI];
1710  label index = oldPointToIndex[oldPointI];
1711  if (index != -1)
1712  {
1713  pointStat[pointI] = classifyFeaturePoint(pointI);
1714  }
1715  }
1716 
1717  // Reset based on new point and edge status
1718  labelList sortedToOriginalPoint;
1719  labelList sortedToOriginalEdge;
1720  setFromStatus
1721  (
1722  pointStat,
1723  edgeStat,
1724  sortedToOriginalPoint,
1725  sortedToOriginalEdge
1726  );
1727 
1728  // Update the overall pointMap, edgeMap
1729  pointMap = labelUIndList(pointMap, sortedToOriginalPoint)();
1730  edgeMap = labelUIndList(edgeMap, sortedToOriginalEdge)();
1731 }
1732 
1733 
1736  const List<extendedEdgeMesh::pointStatus>& pointStat,
1737  const List<extendedEdgeMesh::edgeStatus>& edgeStat,
1738  labelList& sortedToOriginalPoint,
1739  labelList& sortedToOriginalEdge
1740 )
1741 {
1742  // Use pointStatus and edgeStatus to determine new ordering
1743  label pointConcaveStart;
1744  label pointMixedStart;
1745  label pointNonFeatStart;
1746 
1747  label edgeInternalStart;
1748  label edgeFlatStart;
1749  label edgeOpenStart;
1750  label edgeMultipleStart;
1751  sortedOrder
1752  (
1753  pointStat,
1754  edgeStat,
1755  sortedToOriginalPoint,
1756  sortedToOriginalEdge,
1757 
1758  pointConcaveStart,
1759  pointMixedStart,
1760  pointNonFeatStart,
1761 
1762  edgeInternalStart,
1763  edgeFlatStart,
1764  edgeOpenStart,
1765  edgeMultipleStart
1766  );
1767 
1768 
1769  // Order points and edges
1770  labelList reversePointMap(points().size(), -1);
1771  forAll(sortedToOriginalPoint, sortedI)
1772  {
1773  reversePointMap[sortedToOriginalPoint[sortedI]] = sortedI;
1774  }
1775 
1776  edgeList sortedEdges(UIndirectList<edge>(edges(), sortedToOriginalEdge)());
1777  forAll(sortedEdges, sortedI)
1778  {
1779  inplaceRenumber(reversePointMap, sortedEdges[sortedI]);
1780  }
1781 
1782  // Update local data
1783  autoMap
1784  (
1785  pointField(points(), sortedToOriginalPoint),
1786  sortedEdges,
1787  sortedToOriginalPoint,
1788  sortedToOriginalEdge
1789  );
1790 
1791  // Reset the slice starts
1792  concaveStart_ = pointConcaveStart;
1793  mixedStart_ = pointMixedStart;
1794  nonFeatureStart_ = pointNonFeatStart;
1795  internalStart_ = edgeInternalStart;
1796  flatStart_ = edgeFlatStart;
1797  openStart_ = edgeOpenStart;
1798  multipleStart_ = edgeMultipleStart;
1799 }
1800 
1801 
1804  const scalar mergeDist,
1805  labelList& pointMap,
1806  labelList& edgeMap
1807 )
1808 {
1809  const label nOldPoints = points().size();
1810 
1811  // Detect and merge collocated feature points
1812  labelList oldToMerged;
1813  label nNewPoints = ::Foam::mergePoints
1814  (
1815  points(),
1816  SMALL,
1817  false,
1818  oldToMerged
1819  );
1820 
1821  pointMap.setSize(nNewPoints);
1822  pointMap = -1;
1823  forAll(oldToMerged, oldI)
1824  {
1825  label newI = oldToMerged[oldI];
1826  if (pointMap[newI] == -1)
1827  {
1828  pointMap[newI] = oldI;
1829  }
1830  }
1831 
1832  // Renumber edges
1833  edgeList newEdges(edges().size());
1834  forAll(edges(), edgeI)
1835  {
1836  const edge& oldE = edges()[edgeI];
1837  newEdges[edgeI] = edge(oldToMerged[oldE[0]], oldToMerged[oldE[1]]);
1838  }
1839 
1840  // Shuffle basic information (reorders point data)
1841  autoMap
1842  (
1843  pointField(points(), pointMap),
1844  newEdges,
1845  pointMap,
1846  identity(newEdges.size())
1847  );
1848 
1849  // Re-classify the merged points
1850  List<edgeStatus> edgeStat(edges().size());
1851  forAll(edgeStat, edgeI)
1852  {
1853  edgeStat[edgeI] = getEdgeStatus(edgeI);
1854  }
1855 
1856  List<pointStatus> pointStat(points().size());
1857  forAll(pointStat, pointI)
1858  {
1859  pointStat[pointI] = getPointStatus(pointI);
1860  }
1861 
1862  // Re-classify merged points
1863  labelList nPoints(nNewPoints, Zero);
1864  forAll(oldToMerged, oldPointI)
1865  {
1866  nPoints[oldToMerged[oldPointI]]++;
1867  }
1868 
1869  forAll(nPoints, pointI)
1870  {
1871  if (nPoints[pointI] != 1)
1872  {
1873  pointStat[pointI] = classifyFeaturePoint(pointI);
1874  }
1875  }
1876 
1877  labelList sortedToOriginalPoint;
1878  setFromStatus
1879  (
1880  pointStat,
1881  edgeStat,
1882  sortedToOriginalPoint,
1883  edgeMap // point merging above did not affect edge order
1884  );
1885  pointMap = labelUIndList(pointMap, sortedToOriginalPoint)();
1886 
1887  return nNewPoints != nOldPoints;
1888 }
1889 
1890 
1892 {
1893  Info<< nl << "Writing extendedEdgeMesh components to " << prefix
1894  << endl;
1895 
1896  edgeMesh::write(prefix + "_edgeMesh.obj");
1897 
1898  {
1899  OBJstream convexFtPtStr(prefix + "_convexFeaturePts.obj");
1900  Info<< "Writing " << concaveStart_
1901  << " convex feature points to " << convexFtPtStr.name() << endl;
1902 
1903  for(label i = 0; i < concaveStart_; i++)
1904  {
1905  convexFtPtStr.write(points()[i]);
1906  }
1907  }
1908 
1909  {
1910  OBJstream concaveFtPtStr(prefix + "_concaveFeaturePts.obj");
1911  Info<< "Writing " << mixedStart_-concaveStart_
1912  << " concave feature points to "
1913  << concaveFtPtStr.name() << endl;
1914 
1915  for(label i = concaveStart_; i < mixedStart_; i++)
1916  {
1917  concaveFtPtStr.write(points()[i]);
1918  }
1919  }
1920 
1921  {
1922  OBJstream mixedFtPtStr(prefix + "_mixedFeaturePts.obj");
1923  Info<< "Writing " << nonFeatureStart_-mixedStart_
1924  << " mixed feature points to " << mixedFtPtStr.name() << endl;
1925 
1926  for(label i = mixedStart_; i < nonFeatureStart_; i++)
1927  {
1928  mixedFtPtStr.write(points()[i]);
1929  }
1930  }
1931 
1932  {
1933  OBJstream mixedFtPtStructureStr(prefix+"_mixedFeaturePtsStructure.obj");
1934  Info<< "Writing "
1935  << nonFeatureStart_-mixedStart_
1936  << " mixed feature point structure to "
1937  << mixedFtPtStructureStr.name() << endl;
1938 
1939  for(label i = mixedStart_; i < nonFeatureStart_; i++)
1940  {
1941  const labelList& ptEds = pointEdges()[i];
1942 
1943  forAll(ptEds, j)
1944  {
1945  const edge& e = edges()[ptEds[j]];
1946  mixedFtPtStructureStr.write
1947  (
1948  linePointRef(points()[e[0]],
1949  points()[e[1]])
1950  );
1951  }
1952  }
1953  }
1954 
1955  {
1956  OBJstream externalStr(prefix + "_externalEdges.obj");
1957  Info<< "Writing " << internalStart_-externalStart_
1958  << " external edges to " << externalStr.name() << endl;
1959 
1960  for (label i = externalStart_; i < internalStart_; i++)
1961  {
1962  const edge& e = edges()[i];
1963  externalStr.write(linePointRef(points()[e[0]], points()[e[1]]));
1964  }
1965  }
1966 
1967  {
1968  OBJstream internalStr(prefix + "_internalEdges.obj");
1969  Info<< "Writing " << flatStart_-internalStart_
1970  << " internal edges to " << internalStr.name() << endl;
1971 
1972  for (label i = internalStart_; i < flatStart_; i++)
1973  {
1974  const edge& e = edges()[i];
1975  internalStr.write(linePointRef(points()[e[0]], points()[e[1]]));
1976  }
1977  }
1978 
1979  {
1980  OBJstream flatStr(prefix + "_flatEdges.obj");
1981  Info<< "Writing " << openStart_-flatStart_
1982  << " flat edges to " << flatStr.name() << endl;
1983 
1984  for (label i = flatStart_; i < openStart_; i++)
1985  {
1986  const edge& e = edges()[i];
1987  flatStr.write(linePointRef(points()[e[0]], points()[e[1]]));
1988  }
1989  }
1990 
1991  {
1992  OBJstream openStr(prefix + "_openEdges.obj");
1993  Info<< "Writing " << multipleStart_-openStart_
1994  << " open edges to " << openStr.name() << endl;
1995 
1996  for (label i = openStart_; i < multipleStart_; i++)
1997  {
1998  const edge& e = edges()[i];
1999  openStr.write(linePointRef(points()[e[0]], points()[e[1]]));
2000  }
2001  }
2002 
2003  {
2004  OBJstream multipleStr(prefix + "_multipleEdges.obj");
2005  Info<< "Writing " << edges().size()-multipleStart_
2006  << " multiple edges to " << multipleStr.name() << endl;
2007 
2008  for (label i = multipleStart_; i < edges().size(); i++)
2009  {
2010  const edge& e = edges()[i];
2011  multipleStr.write(linePointRef(points()[e[0]], points()[e[1]]));
2012  }
2013  }
2014 
2015  {
2016  OBJstream regionStr(prefix + "_regionEdges.obj");
2017  Info<< "Writing " << regionEdges_.size()
2018  << " region edges to " << regionStr.name() << endl;
2019 
2020  forAll(regionEdges_, i)
2021  {
2022  const edge& e = edges()[regionEdges_[i]];
2023  regionStr.write(linePointRef(points()[e[0]], points()[e[1]]));
2024  }
2025  }
2026 
2027  {
2028  OBJstream edgeDirsStr(prefix + "_edgeDirections.obj");
2029  Info<< "Writing " << edgeDirections_.size()
2030  << " edge directions to " << edgeDirsStr.name() << endl;
2031 
2032  forAll(edgeDirections_, i)
2033  {
2034  const vector& eVec = edgeDirections_[i];
2035  const edge& e = edges()[i];
2036 
2037  edgeDirsStr.write
2038  (
2039  linePointRef(points()[e.start()], eVec + points()[e.start()])
2040  );
2041  }
2042  }
2043 }
2044 
2045 
2047 {
2049 
2050  os << indent << "point classification :" << nl;
2051  os << incrIndent;
2052  os << indent << "convex feature points : "
2053  << setw(8) << concaveStart_-convexStart_
2054  //<< setw(8) << convexStart_
2055  << nl;
2056  os << indent << "concave feature points : "
2057  << setw(8) << mixedStart_-concaveStart_
2058  //<< setw(8) << concaveStart_
2059  << nl;
2060  os << indent << "mixed feature points : "
2061  << setw(8) << nonFeatureStart_-mixedStart_
2062  //<< setw(8) << mixedStart_
2063  << nl;
2064  os << indent << "other (non-feature) points : "
2065  << setw(8) << points().size()-nonFeatureStart_
2066  //<< setw(8) << nonFeatureStart_
2067  << nl;
2068  os << decrIndent;
2069 
2070  os << indent << "edge classification :" << nl;
2071  os << incrIndent;
2072  os << indent << "external (convex angle) edges : "
2073  << setw(8) << internalStart_-externalStart_
2074  //<< setw(8) << externalStart_
2075  << nl;
2076  os << indent << "internal (concave angle) edges : "
2077  << setw(8) << flatStart_-internalStart_
2078  //<< setw(8) << internalStart_
2079  << nl;
2080  os << indent << "flat region edges : "
2081  << setw(8) << openStart_-flatStart_
2082  //<< setw(8) << flatStart_
2083  << nl;
2084  os << indent << "open edges : "
2085  << setw(8) << multipleStart_-openStart_
2086  //<< setw(8) << openStart_
2087  << nl;
2088  os << indent << "multiply connected edges : "
2089  << setw(8) << edges().size()-multipleStart_
2090  //<< setw(8) << multipleStart_
2091  << nl;
2092  os << decrIndent;
2093 }
2094 
2095 
2099  const List<vector>& norms,
2100  const labelList& edNorms,
2101  const vector& fC0tofC1
2102 )
2103 {
2104  label nEdNorms = edNorms.size();
2105 
2106  if (nEdNorms == 1)
2107  {
2108  return OPEN;
2109  }
2110  else if (nEdNorms == 2)
2111  {
2112  const vector& n0(norms[edNorms[0]]);
2113  const vector& n1(norms[edNorms[1]]);
2114 
2115  if ((n0 & n1) > cosNormalAngleTol_)
2116  {
2117  return FLAT;
2118  }
2119  else if ((fC0tofC1 & n0) > 0.0)
2120  {
2121  return INTERNAL;
2122  }
2123  else
2124  {
2125  return EXTERNAL;
2126  }
2127  }
2128  else if (nEdNorms > 2)
2129  {
2130  return MULTIPLE;
2131  }
2132 
2133  // There is a problem - the edge has no normals
2134  return NONE;
2135 }
2136 
2137 
2140  const List<extendedEdgeMesh::pointStatus>& pointStat,
2141  const List<extendedEdgeMesh::edgeStatus>& edgeStat,
2142  labelList& sortedToOriginalPoint,
2143  labelList& sortedToOriginalEdge,
2144 
2145  label& pointConcaveStart,
2146  label& pointMixedStart,
2147  label& pointNonFeatStart,
2148 
2149  label& edgeInternalStart,
2150  label& edgeFlatStart,
2151  label& edgeOpenStart,
2152  label& edgeMultipleStart
2153 )
2154 {
2155  sortedToOriginalPoint.setSize(pointStat.size());
2156  sortedToOriginalPoint = -1;
2157 
2158  sortedToOriginalEdge.setSize(edgeStat.size());
2159  sortedToOriginalEdge = -1;
2160 
2161 
2162  // Order edges
2163  // ~~~~~~~~~~~
2164 
2165  label nConvex = 0;
2166  label nConcave = 0;
2167  label nMixed = 0;
2168  label nNonFeat = 0;
2169 
2170  forAll(pointStat, pointI)
2171  {
2172  switch (pointStat[pointI])
2173  {
2175  nConvex++;
2176  break;
2177 
2179  nConcave++;
2180  break;
2181 
2183  nMixed++;
2184  break;
2185 
2187  nNonFeat++;
2188  break;
2189 
2190  default:
2192  << "Problem" << exit(FatalError);
2193  break;
2194  }
2195  }
2196 
2197  label convexStart = 0;
2198  label concaveStart = nConvex;
2199  label mixedStart = concaveStart+nConcave;
2200  label nonFeatStart = mixedStart+nMixed;
2201 
2202 
2203  // Copy to parameters
2204  pointConcaveStart = concaveStart;
2205  pointMixedStart = mixedStart;
2206  pointNonFeatStart = nonFeatStart;
2207 
2208  forAll(pointStat, pointI)
2209  {
2210  switch (pointStat[pointI])
2211  {
2213  sortedToOriginalPoint[convexStart++] = pointI;
2214  break;
2215 
2217  sortedToOriginalPoint[concaveStart++] = pointI;
2218  break;
2219 
2221  sortedToOriginalPoint[mixedStart++] = pointI;
2222  break;
2223 
2225  sortedToOriginalPoint[nonFeatStart++] = pointI;
2226  break;
2227  }
2228  }
2229 
2230 
2231  // Order edges
2232  // ~~~~~~~~~~~
2233 
2234  label nExternal = 0;
2235  label nInternal = 0;
2236  label nFlat = 0;
2237  label nOpen = 0;
2238  label nMultiple = 0;
2239 
2240  forAll(edgeStat, edgeI)
2241  {
2242  switch (edgeStat[edgeI])
2243  {
2245  nExternal++;
2246  break;
2247 
2249  nInternal++;
2250  break;
2251 
2253  nFlat++;
2254  break;
2255 
2257  nOpen++;
2258  break;
2259 
2261  nMultiple++;
2262  break;
2263 
2265  default:
2267  << "Problem" << exit(FatalError);
2268  break;
2269  }
2270  }
2271 
2272  label externalStart = 0;
2273  label internalStart = nExternal;
2274  label flatStart = internalStart + nInternal;
2275  label openStart = flatStart + nFlat;
2276  label multipleStart = openStart + nOpen;
2277 
2278 
2279  // Copy to parameters
2280  edgeInternalStart = internalStart;
2281  edgeFlatStart = flatStart;
2282  edgeOpenStart = openStart;
2283  edgeMultipleStart = multipleStart;
2284 
2285  forAll(edgeStat, edgeI)
2286  {
2287  switch (edgeStat[edgeI])
2288  {
2290  sortedToOriginalEdge[externalStart++] = edgeI;
2291  break;
2292 
2294  sortedToOriginalEdge[internalStart++] = edgeI;
2295  break;
2296 
2298  sortedToOriginalEdge[flatStart++] = edgeI;
2299  break;
2300 
2302  sortedToOriginalEdge[openStart++] = edgeI;
2303  break;
2304 
2306  sortedToOriginalEdge[multipleStart++] = edgeI;
2307  break;
2308 
2310  default:
2312  << "Problem" << exit(FatalError);
2313  break;
2314  }
2315  }
2316 }
2317 
2318 
2319 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
2320 
2321 Foam::Istream& Foam::operator>>
2323  Istream& is,
2325 )
2326 {
2327  label type;
2328  is >> type;
2329 
2330  vt = static_cast<Foam::extendedEdgeMesh::sideVolumeType>(type);
2331 
2332  is.check(FUNCTION_NAME);
2333  return is;
2334 }
2335 
2336 
2337 Foam::Ostream& Foam::operator<<
2339  Ostream& os,
2341 )
2342 {
2343  os << static_cast<label>(vt);
2344 
2345  return os;
2346 }
2347 
2348 
2350 {
2351  //fileFormats::extendedEdgeMeshFormat::write(os, em.points_, em.edges_);
2352  os << "// points" << nl
2353  << em.points() << nl
2354  << "// edges" << nl
2355  << em.edges() << nl
2356  << "// concaveStart mixedStart nonFeatureStart" << nl
2357  << em.concaveStart_ << token::SPACE
2358  << em.mixedStart_ << token::SPACE
2359  << em.nonFeatureStart_ << nl
2360  << "// internalStart flatStart openStart multipleStart" << nl
2361  << em.internalStart_ << token::SPACE
2362  << em.flatStart_ << token::SPACE
2363  << em.openStart_ << token::SPACE
2364  << em.multipleStart_ << nl
2365  << "// normals" << nl
2366  << em.normals_ << nl
2367  << "// normal volume types" << nl
2368  << em.normalVolumeTypes_ << nl
2369  << "// normalDirections" << nl
2370  << em.normalDirections_ << nl
2371  << "// edgeNormals" << nl
2372  << em.edgeNormals_ << nl
2373  << "// featurePointNormals" << nl
2374  << em.featurePointNormals_ << nl
2375  << "// featurePointEdges" << nl
2376  << em.featurePointEdges_ << nl
2377  << "// regionEdges" << nl
2378  << em.regionEdges_
2379  << endl;
2380 
2382  return os;
2383 }
2384 
2385 
2387 {
2388  //fileFormats::extendedEdgeMeshFormat::read(is, em.points_, em.edges_);
2389  is >> static_cast<edgeMesh&>(em)
2390  >> em.concaveStart_
2391  >> em.mixedStart_
2392  >> em.nonFeatureStart_
2393  >> em.internalStart_
2394  >> em.flatStart_
2395  >> em.openStart_
2396  >> em.multipleStart_
2397  >> em.normals_
2398  >> em.normalVolumeTypes_
2399  >> em.normalDirections_
2400  >> em.edgeNormals_
2401  >> em.featurePointNormals_
2402  >> em.featurePointEdges_
2403  >> em.regionEdges_;
2404 
2405  is.check(FUNCTION_NAME);
2406  return is;
2407 }
2408 
2409 
2410 // ************************************************************************* //
Foam::extendedEdgeMesh::NONE
Unclassified (consistency with surfaceFeatures)
Definition: extendedEdgeMesh.H:111
Foam::word::lessExt
word lessExt() const
Return word without extension (part before last .)
Definition: word.C:113
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::extendedEdgeMesh::regionEdges_
labelList regionEdges_
Feature edges which are on the boundary between regions.
Definition: extendedEdgeMesh.H:192
Foam::edgeMesh::points
const pointField & points() const noexcept
Return points.
Definition: edgeMeshI.H:99
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::extendedEdgeMesh::edgeStatus
edgeStatus
Definition: extendedEdgeMesh.H:104
Foam::searchableSurface::findLineAll
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &) const =0
Get all intersections in order from start to end.
Foam::extendedEdgeMesh::sortedOrder
static void sortedOrder(const List< extendedEdgeMesh::pointStatus > &pointStat, const List< extendedEdgeMesh::edgeStatus > &edgeStat, labelList &sortedToOriginalPoint, labelList &sortedToOriginalEdge, label &pointConcaveStart, label &pointMixedStart, label &pointNonFeatStart, label &edgeInternalStart, label &edgeFlatStart, label &edgeOpenStart, label &edgeMultipleStart)
Determine the ordering.
Definition: extendedEdgeMesh.C:2139
Foam::extendedEdgeMesh::featurePointNormals
const labelListList & featurePointNormals() const
Return the indices of the normals that are adjacent to the.
Definition: extendedEdgeMeshI.H:177
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatch.C:262
Foam::labelMax
constexpr label labelMax
Definition: label.H:61
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::uniform
static Vector< Cmpt > uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:164
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::OBJstream
OFstream that keeps track of vertices.
Definition: OBJstream.H:58
searchableSurface.H
Foam::extendedEdgeMesh::allNearestFeatureEdges
void allNearestFeatureEdges(const point &sample, const scalar searchRadiusSqr, List< pointIndexHit > &info) const
Find all the feature edges within searchDistSqr of sample.
Definition: extendedEdgeMesh.C:771
Foam::extendedEdgeMesh::canReadType
static bool canReadType(const word &fileType, bool verbose=false)
Can we read this file format?
Definition: extendedEdgeMesh.C:113
Foam::OFstream::name
virtual const fileName & name() const
Read/write access to the name of the stream.
Definition: OSstream.H:107
Foam::PointHit
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:53
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< label >
Foam::extendedEdgeMesh::FLAT
Neither concave or convex, on a flat surface.
Definition: extendedEdgeMesh.H:108
Foam::extendedEdgeMesh::edgeTreesByType
const PtrList< indexedOctree< treeDataEdge > > & edgeTreesByType() const
Demand driven construction of octree for boundary edges by type.
Definition: extendedEdgeMesh.C:909
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:54
Foam::extendedEdgeMesh::concaveStart_
label concaveStart_
Index of the start of the concave feature points.
Definition: extendedEdgeMesh.H:145
surfaceFeatures.H
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::extendedEdgeMesh::nearestFeatureEdgeByType
void nearestFeatureEdgeByType(const point &sample, const scalarField &searchDistSqr, List< pointIndexHit > &info) const
Find the nearest point on each type of feature edge.
Definition: extendedEdgeMesh.C:704
Foam::extendedEdgeMesh::mergePointsAndSort
bool mergePointsAndSort(const scalar mergeDist, labelList &pointMap, labelList &edgeMap)
Geometric merge points. Returns true if any points merged.
Definition: extendedEdgeMesh.C:1803
Foam::DynamicField::capacity
label capacity() const noexcept
Size of the underlying storage.
Definition: DynamicFieldI.H:327
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::extendedEdgeMesh::autoMap
void autoMap(const pointField &subPoints, const edgeList &subEdges, const labelList &pointMap, const labelList &edgeMap)
Update with derived geometry.
Definition: extendedEdgeMesh.C:1412
Foam::extendedEdgeMesh::regionEdges
const labelList & regionEdges() const
Return the feature edges which are on the boundary between.
Definition: extendedEdgeMeshI.H:218
Foam::extendedEdgeMesh::normals
const vectorField & normals() const
Return the normals of the surfaces adjacent to the feature edges.
Definition: extendedEdgeMeshI.H:90
Foam::bitSet::set
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:574
Foam::extendedEdgeMesh::edgeTree
const indexedOctree< treeDataEdge > & edgeTree() const
Demand driven construction of octree for boundary edges.
Definition: extendedEdgeMesh.C:867
Foam::objectRegistry::clear
void clear()
Clear all entries from the registry.
Definition: objectRegistry.C:318
Foam::extendedEdgeMesh::NONFEATURE
Not a feature point.
Definition: extendedEdgeMesh.H:99
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:230
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::extendedEdgeMesh::multipleStart
label multipleStart() const
Return the index of the start of the multiply-connected feature.
Definition: extendedEdgeMeshI.H:78
Foam::extendedEdgeMesh::mixedStart_
label mixedStart_
Index of the start of the mixed type feature points.
Definition: extendedEdgeMesh.H:148
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
Foam::DynamicList::capacity
label capacity() const noexcept
Size of the underlying storage.
Definition: DynamicListI.H:287
Foam::word::ext
word ext() const
Return file name extension (part after last .)
Definition: word.C:126
triSurface.H
Foam::bitSet::test
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Definition: bitSetI.H:520
Foam::HashSet< word, Hash< word > >
Foam::OBJstream::write
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:346
Foam::extendedEdgeMesh::normalVolumeTypes_
List< sideVolumeType > normalVolumeTypes_
Type per normal: which side of normal to mesh.
Definition: extendedEdgeMesh.H:170
Foam::edgeMesh::write
static void write(const fileName &name, const edgeMesh &mesh, IOstreamOption streamOpt=IOstreamOption(), const dictionary &options=dictionary::null)
Write to file (format implicit in the extension)
Definition: edgeMeshIO.C:114
Foam::inplaceRenumber
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values (not the indices) of a list.
Definition: ListOpsTemplates.C:61
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
Foam::treeDataPoint
Holds (reference to) pointField. Encapsulation of data needed for octree searches....
Definition: treeDataPoint.H:62
Foam::DynamicList::setCapacity
void setCapacity(const label len)
Alter the size of the underlying storage.
Definition: DynamicListI.H:303
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::extendedEdgeMesh::sideVolumeTypeNames_
static const Enum< sideVolumeType > sideVolumeTypeNames_
Definition: extendedEdgeMesh.H:125
Foam::extendedEdgeMesh::INTERNAL
"Concave" edge
Definition: extendedEdgeMesh.H:107
Foam::DynamicField
Dynamically sized Field.
Definition: DynamicField.H:49
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::extendedEdgeMesh::edgeDirections
const vectorField & edgeDirections() const
Return the edgeDirection vectors.
Definition: extendedEdgeMeshI.H:103
Foam::extendedEdgeMesh::classifyEdge
static edgeStatus classifyEdge(const List< vector > &norms, const labelList &edNorms, const vector &fC0tofC1)
Classify the type of feature edge. Requires face centre 0 to face.
Definition: extendedEdgeMesh.C:2098
Foam::surfaceFeatures::featurePoints
const labelList & featurePoints() const
Return feature point list.
Definition: surfaceFeatures.H:250
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::extendedEdgeMesh::canWriteType
static bool canWriteType(const word &fileType, bool verbose=false)
Can we write this file format type?
Definition: extendedEdgeMesh.C:125
Foam::extendedEdgeMesh::nonFeatureStart
label nonFeatureStart() const
Return the index of the start of the non-feature points.
Definition: extendedEdgeMeshI.H:48
Foam::extendedEdgeMesh::concaveStart
label concaveStart() const
Return the index of the start of the concave feature points.
Definition: extendedEdgeMeshI.H:36
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:52
Foam::extendedEdgeMesh::readTypes
static wordHashSet readTypes()
Summary of supported read file types.
Definition: extendedEdgeMesh.C:101
Foam::DynamicField::setCapacity
void setCapacity(const label len)
Alter the size of the underlying storage.
Definition: DynamicFieldI.H:343
Foam::volumeType
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition: volumeType.H:60
Foam::Field< vector >
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:76
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::extendedEdgeMesh::MIXED
A point surrounded by both convex and concave edges.
Definition: extendedEdgeMesh.H:98
Foam::extendedEdgeMesh::setFromStatus
void setFromStatus(const List< extendedEdgeMesh::pointStatus > &pointStat, const List< extendedEdgeMesh::edgeStatus > &edgeStat, labelList &sortedToOriginalPoint, labelList &sortedToOriginalEdge)
Order according to point and edge status.
Definition: extendedEdgeMesh.C:1735
Foam::extendedEdgeMesh::sideVolumeType
sideVolumeType
Normals point to the outside.
Definition: extendedEdgeMesh.H:117
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
Foam::PointHit::rawPoint
const point_type & rawPoint() const noexcept
The point, no checks.
Definition: PointHit.H:172
Foam::searchableSurface
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Definition: searchableSurface.H:69
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
Foam::surfaceFeatures
Holds feature edges/points of surface.
Definition: surfaceFeatures.H:73
Foam::edgeMesh::writeStats
virtual void writeStats(Ostream &) const
Definition: edgeMeshIO.C:125
Foam::indexedOctree< Foam::treeDataPoint >
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::treeDataEdge
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:56
samples
scalarField samples(nIntervals, Zero)
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::extendedEdgeMesh::multipleStart_
label multipleStart_
Index of the start of the multiply-connected feature edges.
Definition: extendedEdgeMesh.H:163
Foam::extendedEdgeMesh::clear
virtual void clear()
Clear all storage.
Definition: extendedEdgeMesh.C:1006
Foam::extendedEdgeMesh::cosNormalAngleTol_
static scalar cosNormalAngleTol_
Angular closeness tolerance for treating normals as the same.
Definition: extendedEdgeMesh.H:128
extendedEdgeMesh.H
Foam::surfaceFeatures::featureEdges
const labelList & featureEdges() const
Return feature edge list.
Definition: surfaceFeatures.H:256
Foam::extendedEdgeMesh::nearestFeaturePoint
void nearestFeaturePoint(const point &sample, scalar searchDistSqr, pointIndexHit &info) const
Find nearest surface edge for the sample point.
Definition: extendedEdgeMesh.C:653
Foam::extendedEdgeMesh::select
void select(const searchableSurface &surf, const volumeType volType, labelList &pMap, labelList &eMap)
Remove outside/inside edges. volType denotes which side to keep.
Definition: extendedEdgeMesh.C:298
cut
Patchify triangles based on orientation w.r.t other (triangulated or triangulatable) surfaces.
Foam::extendedEdgeMesh::pointStatusNames_
static const Enum< pointStatus > pointStatusNames_
Definition: extendedEdgeMesh.H:102
newPointi
label newPointi
Definition: readKivaGrid.H:496
Foam::extendedEdgeMesh::externalStart_
static label externalStart_
Index of the start of the external feature edges - static as 0.
Definition: extendedEdgeMesh.H:139
Foam::Field::rmap
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:466
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
Foam::FatalError
error FatalError
Foam::extendedEdgeMesh::internalStart
label internalStart() const
Return the index of the start of the internal feature edges.
Definition: extendedEdgeMeshI.H:60
Foam::extendedEdgeMesh::edgeNormals
const labelListList & edgeNormals() const
Return the indices of the normals that are adjacent to the.
Definition: extendedEdgeMeshI.H:146
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::extendedEdgeMesh::MULTIPLE
Multiply connected (connected to more than two faces)
Definition: extendedEdgeMesh.H:110
Foam::fileName::ext
word ext() const
Return file name extension (part after last .)
Definition: fileNameI.H:218
Foam::extendedEdgeMesh::edgeNormals_
labelListList edgeNormals_
Indices of the normals that are adjacent to the feature edges.
Definition: extendedEdgeMesh.H:181
Foam::extendedEdgeMesh::CONVEX
Fully convex point (w.r.t normals)
Definition: extendedEdgeMesh.H:96
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::extendedEdgeMesh::pointStatus
pointStatus
Definition: extendedEdgeMesh.H:94
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:353
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:339
Foam::extendedEdgeMesh::mixedStart
label mixedStart() const
Return the index of the start of the mixed type feature points.
Definition: extendedEdgeMeshI.H:42
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::HashTable::transfer
void transfer(HashTable< T, Key, Hash > &rhs)
Transfer contents into this table.
Definition: HashTable.C:671
Foam::mergePoints
label mergePoints(const PointList &points, const scalar mergeTol, const bool verbose, labelList &pointMap, typename PointList::const_reference origin=PointList::value_type::zero)
Sorts and merges points. All points closer than/equal mergeTol get merged.
Foam::extendedEdgeMesh::edgeStatusNames_
static const Enum< edgeStatus > edgeStatusNames_
Definition: extendedEdgeMesh.H:114
Random.H
Foam::edgeMesh::clear
virtual void clear()
Clear all storage.
Definition: edgeMesh.C:116
Foam::extendedEdgeMesh::normals_
vectorField normals_
Normals of the features, to be referred to by index by both feature.
Definition: extendedEdgeMesh.H:167
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::extendedEdgeMesh::transfer
void transfer(extendedEdgeMesh &mesh)
Transfer the contents of the argument and annul the argument.
Definition: extendedEdgeMesh.C:974
Foam::surfaceFeatures::surface
const triSurface & surface() const
Definition: surfaceFeatures.H:244
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::extendedEdgeMesh::internalStart_
label internalStart_
Index of the start of the internal feature edges.
Definition: extendedEdgeMesh.H:154
Foam::extendedEdgeMesh::convexStart_
static label convexStart_
Index of the start of the convex feature points - static as 0.
Definition: extendedEdgeMesh.H:136
Foam::extendedEdgeMesh::normalDirections_
labelListList normalDirections_
Starting directions for the edges.
Definition: extendedEdgeMesh.H:178
Foam::extendedEdgeMesh::openStart
label openStart() const
Return the index of the start of the open feature edges.
Definition: extendedEdgeMeshI.H:72
Time.H
Foam::edgeMesh::edges
const edgeList & edges() const noexcept
Return edges.
Definition: edgeMeshI.H:105
Foam::extendedEdgeMesh
Description of feature edges and points.
Definition: extendedEdgeMesh.H:85
Foam::linePointRef
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:47
Foam::extendedEdgeMesh::pointTree
const indexedOctree< treeDataPoint > & pointTree() const
Demand driven construction of octree for feature points.
Definition: extendedEdgeMesh.C:827
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
DynamicField.H
Foam::extendedEdgeMesh::extendedEdgeMesh
extendedEdgeMesh()
Default construct.
Definition: extendedEdgeMesh.C:377
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::extendedEdgeMesh::OPEN
Only connected to a single face.
Definition: extendedEdgeMesh.H:109
Foam::extendedEdgeMesh::add
void add(const extendedEdgeMesh &fem)
Add extendedEdgeMesh. No filtering of duplicates.
Definition: extendedEdgeMesh.C:1030
Foam::Vector< scalar >
Foam::extendedEdgeMesh::trim
void trim(const searchableSurface &surf, const volumeType volType, labelList &pointMap, labelList &edgeMap)
Trim to surface. Keep volType side. Return map from current back.
Definition: extendedEdgeMesh.C:1652
Foam::List< label >
Foam::extendedEdgeMesh::writeStats
virtual void writeStats(Ostream &os) const
Dump some information.
Definition: extendedEdgeMesh.C:2046
Foam::extendedEdgeMesh::EXTERNAL
"Convex" edge
Definition: extendedEdgeMesh.H:106
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::token::SPACE
Space [isspace].
Definition: token.H:125
Foam::UList< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::extendedEdgeMesh::flipNormals
void flipNormals()
Flip normals. All concave become convex, all internal external.
Definition: extendedEdgeMesh.C:1300
Foam::extendedEdgeMesh::nearestFeatureEdge
void nearestFeatureEdge(const point &sample, scalar searchDistSqr, pointIndexHit &info) const
Find nearest surface edge for the sample point.
Definition: extendedEdgeMesh.C:668
Foam::volumeType::INSIDE
A location inside the volume.
Definition: volumeType.H:68
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:295
Foam::extendedEdgeMesh::classifyFeaturePoint
pointStatus classifyFeaturePoint(label ptI) const
Classify the type of feature point. Requires valid stored member.
Definition: extendedEdgeMesh.C:152
Foam::extendedEdgeMesh::cut
void cut(const searchableSurface &, labelList &pMap, labelList &eMap, labelList &pointsFromEdge, labelList &oldEdge, labelList &surfTri)
Cut edges with surface. Return map from cut points&edges back.
Definition: extendedEdgeMesh.C:196
Foam::surfaceFeatures::nRegionEdges
label nRegionEdges() const
Return number of region edges.
Definition: surfaceFeatures.H:274
Foam::wordHashSet
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:77
Foam::extendedEdgeMesh::featurePointEdges_
labelListList featurePointEdges_
Indices of feature edges attached to feature points. The edges are.
Definition: extendedEdgeMesh.H:189
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
rndGen
Random rndGen
Definition: createFields.H:23
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::extendedEdgeMesh::writeTypes
static wordHashSet writeTypes()
Summary of supported write file types.
Definition: extendedEdgeMesh.C:107
Foam::extendedEdgeMesh::flatStart_
label flatStart_
Index of the start of the flat feature edges.
Definition: extendedEdgeMesh.H:157
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::sortedOrder
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::one::minus
A class representing the concept of -1.
Definition: one.H:105
Foam::extendedEdgeMesh::nonFeatureStart_
label nonFeatureStart_
Index of the start of the non-feature points.
Definition: extendedEdgeMesh.H:151
Foam::extendedEdgeMesh::read
bool read(const fileName &name, const word &ext)
Read from file. Chooses reader based on explicit extension.
Definition: extendedEdgeMesh.C:641
Foam::extendedEdgeMesh::writeObj
void writeObj(const fileName &prefix) const
Write all components of the extendedEdgeMesh as obj files.
Definition: extendedEdgeMesh.C:1891
Foam::extendedEdgeMesh::allNearestFeaturePoints
void allNearestFeaturePoints(const point &sample, scalar searchRadiusSqr, List< pointIndexHit > &info) const
Find all the feature points within searchDistSqr of sample.
Definition: extendedEdgeMesh.C:740
Foam::DynamicField::append
DynamicField< T, SizeMin > & append(const T &val)
Append an element at the end of the list.
Definition: DynamicFieldI.H:587
Foam::extendedEdgeMesh::openStart_
label openStart_
Index of the start of the open feature edges.
Definition: extendedEdgeMesh.H:160
triSurfaceMesh.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::edgeMesh
Mesh data needed to do the Finite Area discretisation.
Definition: edgeFaMesh.H:53
Foam::PointHit::hit
bool hit() const noexcept
Is there a hit.
Definition: PointHit.H:121
edgeMeshFormatsCore.H
Foam::edgeMesh::transfer
void transfer(edgeMesh &mesh)
Transfer the contents of the argument and annul the argument.
Definition: edgeMesh.C:124
OBJstream.H
Foam::extendedEdgeMesh::flatStart
label flatStart() const
Return the index of the start of the flat feature edges.
Definition: extendedEdgeMeshI.H:66
Foam::extendedEdgeMesh::featurePointNormals_
labelListList featurePointNormals_
Indices of the normals that are adjacent to the feature points.
Definition: extendedEdgeMesh.H:185
Foam::volumeType::OUTSIDE
A location outside the volume.
Definition: volumeType.H:69
Foam::searchableSurface::getVolumeType
virtual void getVolumeType(const pointField &, List< volumeType > &) const =0
Determine type (inside/outside) for point.
Foam::extendedEdgeMesh::canRead
static bool canRead(const fileName &name, bool verbose=false)
Can we read this file format?
Definition: extendedEdgeMesh.C:137
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79
sample
Minimal example by using system/controlDict.functions:
Foam::extendedEdgeMesh::CONCAVE
Fully concave point.
Definition: extendedEdgeMesh.H:97
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:58
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265