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