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