extendedEdgeMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "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 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
41 
42  template<>
43  const char* Foam::NamedEnum
44  <
46  4
47  >::names[] =
48  {
49  "convex",
50  "concave",
51  "mixed",
52  "nonFeature"
53  };
54 
55  template<>
56  const char* Foam::NamedEnum
57  <
59  6
60  >::names[] =
61  {
62  "external",
63  "internal",
64  "flat",
65  "open",
66  "multiple",
67  "none"
68  };
69 
70  template<>
71  const char* Foam::NamedEnum
72  <
74  4
75  >::names[] =
76  {
77  "inside",
78  "outside",
79  "both",
80  "neither"
81  };
82 }
83 
86 
89 
92 
94  Foam::cos(degToRad(0.1));
95 
96 
98 
99 
101 
102 
104 
105 
107 
108 
110 {
111  return wordHashSet(*fileExtensionConstructorTablePtr_);
112 }
113 
114 
116 {
117  return wordHashSet(*writefileExtensionMemberFunctionTablePtr_);
118 }
119 
120 
121 
122 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
123 
125 (
126  const word& ext,
127  const bool verbose
128 )
129 {
130  return edgeMeshFormatsCore::checkSupport
131  (
132  readTypes(),
133  ext,
134  verbose,
135  "reading"
136  );
137 }
138 
139 
141 (
142  const word& ext,
143  const bool verbose
144 )
145 {
146  return edgeMeshFormatsCore::checkSupport
147  (
148  writeTypes(),
149  ext,
150  verbose,
151  "writing"
152  );
153 }
154 
155 
157 (
158  const fileName& name,
159  const bool verbose
160 )
161 {
162  word ext = name.ext();
163  if (ext == "gz")
164  {
165  ext = name.lessExt().ext();
166  }
167  return canReadType(ext, verbose);
168 }
169 
170 
171 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
172 
175 (
176  label ptI
177 ) const
178 {
179  const labelList& ptEds(pointEdges()[ptI]);
180 
181  label nPtEds = ptEds.size();
182  label nExternal = 0;
183  label nInternal = 0;
184 
185  if (nPtEds == 0)
186  {
187  // There are no edges attached to the point, this is a problem
188  return NONFEATURE;
189  }
190 
191  forAll(ptEds, i)
192  {
193  edgeStatus edStat = getEdgeStatus(ptEds[i]);
194 
195  if (edStat == EXTERNAL)
196  {
197  nExternal++;
198  }
199  else if (edStat == INTERNAL)
200  {
201  nInternal++;
202  }
203  }
204 
205  if (nExternal == nPtEds)
206  {
207  return CONVEX;
208  }
209  else if (nInternal == nPtEds)
210  {
211  return CONCAVE;
212  }
213  else
214  {
215  return MIXED;
216  }
217 }
218 
219 
222 (
223  const List<vector>& norms,
224  const labelList& edNorms,
225  const vector& fC0tofC1
226 )
227 {
228  label nEdNorms = edNorms.size();
229 
230  if (nEdNorms == 1)
231  {
232  return OPEN;
233  }
234  else if (nEdNorms == 2)
235  {
236  const vector& n0(norms[edNorms[0]]);
237  const vector& n1(norms[edNorms[1]]);
238 
239  if ((n0 & n1) > cosNormalAngleTol_)
240  {
241  return FLAT;
242  }
243  else if ((fC0tofC1 & n0) > 0.0)
244  {
245  return INTERNAL;
246  }
247  else
248  {
249  return EXTERNAL;
250  }
251  }
252  else if (nEdNorms > 2)
253  {
254  return MULTIPLE;
255  }
256  else
257  {
258  // There is a problem - the edge has no normals
259  return NONE;
260  }
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
265 
267 :
268  edgeMesh(pointField(0), edgeList(0)),
269  concaveStart_(0),
270  mixedStart_(0),
271  nonFeatureStart_(0),
272  internalStart_(0),
273  flatStart_(0),
274  openStart_(0),
275  multipleStart_(0),
276  normals_(0),
277  normalVolumeTypes_(0),
278  edgeDirections_(0),
279  normalDirections_(0),
280  edgeNormals_(0),
281  featurePointNormals_(0),
282  featurePointEdges_(0),
283  regionEdges_(0),
284  pointTree_(),
285  edgeTree_(),
286  edgeTreesByType_()
287 {}
288 
289 
291 :
292  edgeMesh(fem),
293  concaveStart_(fem.concaveStart()),
294  mixedStart_(fem.mixedStart()),
295  nonFeatureStart_(fem.nonFeatureStart()),
296  internalStart_(fem.internalStart()),
297  flatStart_(fem.flatStart()),
298  openStart_(fem.openStart()),
299  multipleStart_(fem.multipleStart()),
300  normals_(fem.normals()),
301  normalVolumeTypes_(fem.normalVolumeTypes()),
302  edgeDirections_(fem.edgeDirections()),
303  normalDirections_(fem.normalDirections()),
304  edgeNormals_(fem.edgeNormals()),
305  featurePointNormals_(fem.featurePointNormals()),
306  featurePointEdges_(fem.featurePointEdges()),
307  regionEdges_(fem.regionEdges()),
308  pointTree_(),
309  edgeTree_(),
310  edgeTreesByType_()
311 {}
312 
313 
315 :
316  edgeMesh(move(fem)),
317  concaveStart_(fem.concaveStart()),
318  mixedStart_(fem.mixedStart()),
319  nonFeatureStart_(fem.nonFeatureStart()),
320  internalStart_(fem.internalStart()),
321  flatStart_(fem.flatStart()),
322  openStart_(fem.openStart()),
323  multipleStart_(fem.multipleStart()),
324  normals_(move(fem.normals())),
325  normalVolumeTypes_(move(fem.normalVolumeTypes())),
326  edgeDirections_(move(fem.edgeDirections())),
327  normalDirections_(move(fem.normalDirections())),
328  edgeNormals_(move(fem.edgeNormals())),
329  featurePointNormals_(move(fem.featurePointNormals())),
330  featurePointEdges_(move(fem.featurePointEdges())),
331  regionEdges_(move(fem.regionEdges())),
332  pointTree_(),
333  edgeTree_(),
334  edgeTreesByType_()
335 {}
336 
337 
339 {
340  is >> *this;
341 }
342 
343 
345 (
346  pointField&& pointLst,
347  edgeList&& edgeLst
348 )
349 :
350  edgeMesh(move(pointLst), move(edgeLst)),
351  concaveStart_(0),
352  mixedStart_(0),
353  nonFeatureStart_(0),
354  internalStart_(0),
355  flatStart_(0),
356  openStart_(0),
357  multipleStart_(0),
358  normals_(0),
359  normalVolumeTypes_(0),
360  edgeDirections_(0),
361  normalDirections_(0),
362  edgeNormals_(0),
363  featurePointNormals_(0),
364  featurePointEdges_(0),
365  regionEdges_(0),
366  pointTree_(),
367  edgeTree_(),
368  edgeTreesByType_()
369 {}
370 
371 
373 (
374  const surfaceFeatures& sFeat,
375  const boolList& surfBaffleRegions
376 )
377 :
378  edgeMesh(pointField(0), edgeList(0)),
379  concaveStart_(-1),
380  mixedStart_(-1),
381  nonFeatureStart_(-1),
382  internalStart_(-1),
383  flatStart_(-1),
384  openStart_(-1),
385  multipleStart_(-1),
386  normals_(0),
387  normalVolumeTypes_(0),
388  edgeDirections_(0),
389  normalDirections_(0),
390  edgeNormals_(0),
391  featurePointNormals_(0),
392  featurePointEdges_(0),
393  regionEdges_(0),
394  pointTree_(),
395  edgeTree_(),
396  edgeTreesByType_()
397 {
398  // Extract and reorder the data from surfaceFeatures
399  const triSurface& surf = sFeat.surface();
400  const labelList& featureEdges = sFeat.featureEdges();
401  const labelList& featurePoints = sFeat.featurePoints();
402 
403  // Get a labelList of all the featureEdges that are region edges
404  const labelList regionFeatureEdges(identityMap(sFeat.nRegionEdges()));
405 
407  (
408  surf,
409  featureEdges,
410  regionFeatureEdges,
411  featurePoints
412  );
413 
414  const labelListList& edgeFaces = surf.edgeFaces();
415 
417 
418  // Noting when the normal of a face has been used so not to duplicate
419  labelList faceMap(surf.size(), -1);
420 
421  label nAdded = 0;
422 
423  forAll(featureEdges, i)
424  {
425  label sFEI = featureEdges[i];
426 
427  // Pick up the faces adjacent to the feature edge
428  const labelList& eFaces = edgeFaces[sFEI];
429 
430  forAll(eFaces, j)
431  {
432  label eFI = eFaces[j];
433 
434  // Check to see if the points have been already used
435  if (faceMap[eFI] == -1)
436  {
437  normalVolumeTypes_[nAdded++] =
438  (
439  surfBaffleRegions[surf[eFI].region()]
440  ? BOTH
441  : INSIDE
442  );
443 
444  faceMap[eFI] = nAdded - 1;
445  }
446  }
447  }
448 }
449 
450 
452 (
454  const labelList& featureEdges,
455  const labelList& regionFeatureEdges,
456  const labelList& featurePoints
457 )
458 :
459  edgeMesh(pointField(0), edgeList(0)),
460  concaveStart_(-1),
461  mixedStart_(-1),
462  nonFeatureStart_(-1),
463  internalStart_(-1),
464  flatStart_(-1),
465  openStart_(-1),
466  multipleStart_(-1),
467  normals_(0),
468  normalVolumeTypes_(0),
469  edgeDirections_(0),
470  normalDirections_(0),
471  edgeNormals_(0),
472  featurePointNormals_(0),
473  featurePointEdges_(0),
474  regionEdges_(0),
475  pointTree_(),
476  edgeTree_(),
477  edgeTreesByType_()
478 {
480  (
481  surf,
482  featureEdges,
483  regionFeatureEdges,
484  featurePoints
485  );
486 }
487 
488 
490 (
491  const pointField& pts,
492  const edgeList& eds,
493  label concaveStart,
494  label mixedStart,
495  label nonFeatureStart,
496  label internalStart,
497  label flatStart,
498  label openStart,
499  label multipleStart,
500  const vectorField& normals,
501  const List<sideVolumeType>& normalVolumeTypes,
502  const vectorField& edgeDirections,
503  const labelListList& normalDirections,
504  const labelListList& edgeNormals,
505  const labelListList& featurePointNormals,
506  const labelListList& featurePointEdges,
507  const labelList& regionEdges
508 )
509 :
510  edgeMesh(pts, eds),
511  concaveStart_(concaveStart),
512  mixedStart_(mixedStart),
513  nonFeatureStart_(nonFeatureStart),
514  internalStart_(internalStart),
515  flatStart_(flatStart),
516  openStart_(openStart),
517  multipleStart_(multipleStart),
518  normals_(normals),
519  normalVolumeTypes_(normalVolumeTypes),
520  edgeDirections_(edgeDirections),
521  normalDirections_(normalDirections),
522  edgeNormals_(edgeNormals),
523  featurePointNormals_(featurePointNormals),
524  featurePointEdges_(featurePointEdges),
525  regionEdges_(regionEdges),
526  pointTree_(),
527  edgeTree_(),
528  edgeTreesByType_()
529 {}
530 
531 
533 (
534  const fileName& name,
535  const word& ext
536 )
537 :
538  edgeMesh(pointField(0), edgeList(0)),
539  concaveStart_(0),
540  mixedStart_(0),
541  nonFeatureStart_(0),
542  internalStart_(0),
543  flatStart_(0),
544  openStart_(0),
545  multipleStart_(0),
546  normals_(0),
547  normalVolumeTypes_(0),
548  edgeDirections_(0),
549  normalDirections_(0),
550  edgeNormals_(0),
551  featurePointNormals_(0),
552  featurePointEdges_(0),
553  regionEdges_(0),
554  pointTree_(),
555  edgeTree_(),
556  edgeTreesByType_()
557 {
558  read(name, ext);
559 }
560 
561 
563 :
564  edgeMesh(pointField(0), edgeList(0)),
565  concaveStart_(0),
566  mixedStart_(0),
567  nonFeatureStart_(0),
568  internalStart_(0),
569  flatStart_(0),
570  openStart_(0),
571  multipleStart_(0),
572  normals_(0),
573  normalVolumeTypes_(0),
574  edgeDirections_(0),
575  normalDirections_(0),
576  edgeNormals_(0),
577  featurePointNormals_(0),
578  featurePointEdges_(0),
579  regionEdges_(0),
580  pointTree_(),
581  edgeTree_(),
582  edgeTreesByType_()
583 {
584  read(name);
585 }
586 
587 
588 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
589 
591 {}
592 
593 
594 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
595 
597 {
598  word ext = name.ext();
599  if (ext == "gz")
600  {
601  fileName unzipName = name.lessExt();
602  return read(unzipName, unzipName.ext());
603  }
604  else
605  {
606  return read(name, ext);
607  }
608 }
609 
610 
611 // Read from file in given format
613 (
614  const fileName& name,
615  const word& ext
616 )
617 {
618  // read via selector mechanism
619  transfer(New(name, ext)());
620  return true;
621 }
622 
623 
625 (
626  const point& sample,
627  scalar searchDistSqr,
628  pointIndexHit& info
629 ) const
630 {
631  info = pointTree().findNearest
632  (
633  sample,
634  searchDistSqr
635  );
636 }
637 
638 
640 (
641  const point& sample,
642  scalar searchDistSqr,
643  pointIndexHit& info
644 ) const
645 {
646  info = edgeTree().findNearest
647  (
648  sample,
649  searchDistSqr
650  );
651 }
652 
653 
655 (
656  const pointField& samples,
657  const scalarField& searchDistSqr,
658  pointIndexHitList& info
659 ) const
660 {
661  info.setSize(samples.size());
662 
663  forAll(samples, i)
664  {
665  nearestFeatureEdge
666  (
667  samples[i],
668  searchDistSqr[i],
669  info[i]
670  );
671  }
672 }
673 
674 
676 (
677  const point& sample,
678  const scalarField& searchDistSqr,
679  pointIndexHitList& info
680 ) const
681 {
682  const PtrList<indexedOctree<treeDataEdge>>& edgeTrees = edgeTreesByType();
683 
684  info.setSize(edgeTrees.size());
685 
686  labelList sliceStarts(edgeTrees.size());
687 
688  sliceStarts[0] = externalStart_;
689  sliceStarts[1] = internalStart_;
690  sliceStarts[2] = flatStart_;
691  sliceStarts[3] = openStart_;
692  sliceStarts[4] = multipleStart_;
693 
694  forAll(edgeTrees, i)
695  {
696  info[i] = edgeTrees[i].findNearest
697  (
698  sample,
699  searchDistSqr[i]
700  );
701 
702  // The index returned by the indexedOctree is local to the slice of
703  // edges it was supplied with, return the index to the value in the
704  // complete edge list
705 
706  info[i].setIndex(info[i].index() + sliceStarts[i]);
707  }
708 }
709 
710 
712 (
713  const point& sample,
714  scalar searchRadiusSqr,
715  pointIndexHitList& info
716 ) const
717 {
718  // Pick up all the feature points that intersect the search sphere
719  labelList elems = pointTree().findSphere
720  (
721  sample,
722  searchRadiusSqr
723  );
724 
725  DynamicList<pointIndexHit> dynPointHit(elems.size());
726 
727  forAll(elems, elemI)
728  {
729  label index = elems[elemI];
730  label ptI = pointTree().shapes().pointLabels()[index];
731  const point& pt = points()[ptI];
732 
733  pointIndexHit nearHit(true, pt, index);
734 
735  dynPointHit.append(nearHit);
736  }
737 
738  info.transfer(dynPointHit);
739 }
740 
741 
743 (
744  const point& sample,
745  const scalar searchRadiusSqr,
746  pointIndexHitList& info
747 ) const
748 {
749  const PtrList<indexedOctree<treeDataEdge>>& edgeTrees = edgeTreesByType();
750 
751  info.setSize(edgeTrees.size());
752 
753  labelList sliceStarts(edgeTrees.size());
754 
755  sliceStarts[0] = externalStart_;
756  sliceStarts[1] = internalStart_;
757  sliceStarts[2] = flatStart_;
758  sliceStarts[3] = openStart_;
759  sliceStarts[4] = multipleStart_;
760 
761  DynamicList<pointIndexHit> dynEdgeHit(edgeTrees.size()*3);
762 
763  // Loop over all the feature edge types
764  forAll(edgeTrees, i)
765  {
766  // Pick up all the edges that intersect the search sphere
767  labelList elems = edgeTrees[i].findSphere
768  (
769  sample,
770  searchRadiusSqr
771  );
772 
773  forAll(elems, elemI)
774  {
775  label index = elems[elemI];
776  label edgeI = edgeTrees[i].shapes().edgeLabels()[index];
777  const edge& e = edges()[edgeI];
778 
779  pointHit hitPoint = e.line(points()).nearestDist(sample);
780 
781  label hitIndex = index + sliceStarts[i];
782 
783  pointIndexHit nearHit
784  (
785  hitPoint.hit(),
786  hitPoint.rawPoint(),
787  hitIndex
788  );
789 
790  dynEdgeHit.append(nearHit);
791  }
792  }
793 
794  info.transfer(dynEdgeHit);
795 }
796 
797 
799 (
800  const pointIndexHitList& hitList
801 ) const
802 {
803  scalar minDist = GREAT;
804 
805  for
806  (
807  label hi1 = 0;
808  hi1 < hitList.size() - 1;
809  ++hi1
810  )
811  {
812  const pointIndexHit& pHit1 = hitList[hi1];
813 
814  if (pHit1.hit())
815  {
816  const edge& e1 = edges()[pHit1.index()];
817 
818  for
819  (
820  label hi2 = hi1 + 1;
821  hi2 < hitList.size();
822  ++hi2
823  )
824  {
825  const pointIndexHit& pHit2 = hitList[hi2];
826 
827  if (pHit2.hit())
828  {
829  const edge& e2 = edges()[pHit2.index()];
830 
831  // Don't refine if the edges are connected to each other
832  if (!e1.connected(e2))
833  {
834  scalar curDist =
835  mag(pHit1.hitPoint() - pHit2.hitPoint());
836 
837  minDist = min(curDist, minDist);
838  }
839  }
840  }
841  }
842  }
843 
844  return minDist;
845 }
846 
847 
850 {
851  if (pointTree_.empty())
852  {
853  // Slightly extended bb. Slightly off-centred just so on symmetric
854  // geometry there are less face/edge aligned items.
855  treeBoundBox bb(treeBoundBox(points()).extend(1e-4));
856 
857  const labelList featurePointLabels = identityMap(nonFeatureStart_);
858 
859  pointTree_.reset
860  (
862  (
864  (
865  points(),
866  featurePointLabels
867  ),
868  bb, // bb
869  8, // maxLevel
870  10, // leafsize
871  3.0 // duplicity
872  )
873  );
874  }
875 
876  return pointTree_();
877 }
878 
879 
882 {
883  if (edgeTree_.empty())
884  {
885  // Slightly extended bb. Slightly off-centred just so on symmetric
886  // geometry there are less face/edge aligned items.
887  treeBoundBox bb(treeBoundBox(points()).extend(1e-4));
888 
889  labelList allEdges(identityMap(edges().size()));
890 
891  edgeTree_.reset
892  (
894  (
896  (
897  false, // cachebb
898  edges(), // edges
899  points(), // points
900  allEdges // selected edges
901  ),
902  bb, // bb
903  8, // maxLevel
904  10, // leafsize
905  3.0 // duplicity
906  )
907  );
908  }
909 
910  return edgeTree_();
911 }
912 
913 
916 {
917  if (edgeTreesByType_.size() == 0)
918  {
919  edgeTreesByType_.setSize(nEdgeTypes);
920 
921  // Slightly extended bb. Slightly off-centred just so on symmetric
922  // geometry there are less face/edge aligned items.
923  treeBoundBox bb(treeBoundBox(points()).extend(1e-4));
924 
925  labelListList sliceEdges(nEdgeTypes);
926 
927  // External edges
928  sliceEdges[0] =
929  identityMap(internalStart_ - externalStart_) + externalStart_;
930 
931  // Internal edges
932  sliceEdges[1] =
933  identityMap(flatStart_ - internalStart_) + internalStart_;
934 
935  // Flat edges
936  sliceEdges[2] = identityMap(openStart_ - flatStart_) + flatStart_;
937 
938  // Open edges
939  sliceEdges[3] = identityMap(multipleStart_ - openStart_) + openStart_;
940 
941  // Multiple edges
942  sliceEdges[4] =
943  identityMap(edges().size() - multipleStart_) + multipleStart_;
944 
945  forAll(edgeTreesByType_, i)
946  {
947  edgeTreesByType_.set
948  (
949  i,
951  (
953  (
954  false, // cachebb
955  edges(), // edges
956  points(), // points
957  sliceEdges[i] // selected edges
958  ),
959  bb, // bb
960  8, // maxLevel
961  10, // leafsize
962  3.0 // duplicity
963  )
964  );
965  }
966  }
967 
968  return edgeTreesByType_;
969 }
970 
971 
973 {
974  edgeMesh::transfer(mesh);
975 
976  concaveStart_ = mesh.concaveStart_;
977  mixedStart_ = mesh.mixedStart_;
978  nonFeatureStart_ = mesh.nonFeatureStart_;
979  internalStart_ = mesh.internalStart_;
980  flatStart_ = mesh.flatStart_;
981  openStart_ = mesh.openStart_;
982  multipleStart_ = mesh.multipleStart_;
983  normals_.transfer(mesh.normals_);
984  normalVolumeTypes_.transfer(mesh.normalVolumeTypes_);
985  edgeDirections_.transfer(mesh.edgeDirections_);
986  normalDirections_.transfer(mesh.normalDirections_);
987  edgeNormals_.transfer(mesh.edgeNormals_);
988  featurePointNormals_.transfer(mesh.featurePointNormals_);
989  featurePointEdges_.transfer(mesh.featurePointEdges_);
990  regionEdges_.transfer(mesh.regionEdges_);
991  pointTree_ = mesh.pointTree_;
992  edgeTree_ = mesh.edgeTree_;
993  edgeTreesByType_.transfer(mesh.edgeTreesByType_);
994 }
995 
996 
998 {
999  edgeMesh::clear();
1000  concaveStart_ = 0;
1001  mixedStart_ = 0;
1002  nonFeatureStart_ = 0;
1003  internalStart_ = 0;
1004  flatStart_ = 0;
1005  openStart_ = 0;
1006  multipleStart_ = 0;
1007  normals_.clear();
1008  normalVolumeTypes_.clear();
1009  edgeDirections_.clear();
1010  normalDirections_.clear();
1011  edgeNormals_.clear();
1012  featurePointNormals_.clear();
1013  featurePointEdges_.clear();
1014  regionEdges_.clear();
1015  pointTree_.clear();
1016  edgeTree_.clear();
1017  edgeTreesByType_.clear();
1018 }
1019 
1020 
1022 {
1023  // Points
1024  // ~~~~~~
1025 
1026  // From current points into combined points
1027  labelList reversePointMap(points().size());
1028  labelList reverseFemPointMap(fem.points().size());
1029 
1030  label newPointi = 0;
1031  for (label i = 0; i < concaveStart(); i++)
1032  {
1033  reversePointMap[i] = newPointi++;
1034  }
1035  for (label i = 0; i < fem.concaveStart(); i++)
1036  {
1037  reverseFemPointMap[i] = newPointi++;
1038  }
1039 
1040  // Concave
1041  label newConcaveStart = newPointi;
1042  for (label i = concaveStart(); i < mixedStart(); i++)
1043  {
1044  reversePointMap[i] = newPointi++;
1045  }
1046  for (label i = fem.concaveStart(); i < fem.mixedStart(); i++)
1047  {
1048  reverseFemPointMap[i] = newPointi++;
1049  }
1050 
1051  // Mixed
1052  label newMixedStart = newPointi;
1053  for (label i = mixedStart(); i < nonFeatureStart(); i++)
1054  {
1055  reversePointMap[i] = newPointi++;
1056  }
1057  for (label i = fem.mixedStart(); i < fem.nonFeatureStart(); i++)
1058  {
1059  reverseFemPointMap[i] = newPointi++;
1060  }
1061 
1062  // Non-feature
1063  label newNonFeatureStart = newPointi;
1064  for (label i = nonFeatureStart(); i < points().size(); i++)
1065  {
1066  reversePointMap[i] = newPointi++;
1067  }
1068  for (label i = fem.nonFeatureStart(); i < fem.points().size(); i++)
1069  {
1070  reverseFemPointMap[i] = newPointi++;
1071  }
1072 
1073  pointField newPoints(newPointi);
1074  newPoints.rmap(points(), reversePointMap);
1075  newPoints.rmap(fem.points(), reverseFemPointMap);
1076 
1077 
1078  // Edges
1079  // ~~~~~
1080 
1081  // From current edges into combined edges
1082  labelList reverseEdgeMap(edges().size());
1083  labelList reverseFemEdgeMap(fem.edges().size());
1084 
1085  // External
1086  label newEdgeI = 0;
1087  for (label i = 0; i < internalStart(); i++)
1088  {
1089  reverseEdgeMap[i] = newEdgeI++;
1090  }
1091  for (label i = 0; i < fem.internalStart(); i++)
1092  {
1093  reverseFemEdgeMap[i] = newEdgeI++;
1094  }
1095 
1096  // Internal
1097  label newInternalStart = newEdgeI;
1098  for (label i = internalStart(); i < flatStart(); i++)
1099  {
1100  reverseEdgeMap[i] = newEdgeI++;
1101  }
1102  for (label i = fem.internalStart(); i < fem.flatStart(); i++)
1103  {
1104  reverseFemEdgeMap[i] = newEdgeI++;
1105  }
1106 
1107  // Flat
1108  label newFlatStart = newEdgeI;
1109  for (label i = flatStart(); i < openStart(); i++)
1110  {
1111  reverseEdgeMap[i] = newEdgeI++;
1112  }
1113  for (label i = fem.flatStart(); i < fem.openStart(); i++)
1114  {
1115  reverseFemEdgeMap[i] = newEdgeI++;
1116  }
1117 
1118  // Open
1119  label newOpenStart = newEdgeI;
1120  for (label i = openStart(); i < multipleStart(); i++)
1121  {
1122  reverseEdgeMap[i] = newEdgeI++;
1123  }
1124  for (label i = fem.openStart(); i < fem.multipleStart(); i++)
1125  {
1126  reverseFemEdgeMap[i] = newEdgeI++;
1127  }
1128 
1129  // Multiple
1130  label newMultipleStart = newEdgeI;
1131  for (label i = multipleStart(); i < edges().size(); i++)
1132  {
1133  reverseEdgeMap[i] = newEdgeI++;
1134  }
1135  for (label i = fem.multipleStart(); i < fem.edges().size(); i++)
1136  {
1137  reverseFemEdgeMap[i] = newEdgeI++;
1138  }
1139 
1140  edgeList newEdges(newEdgeI);
1141  forAll(edges(), i)
1142  {
1143  const edge& e = edges()[i];
1144  newEdges[reverseEdgeMap[i]] = edge
1145  (
1146  reversePointMap[e[0]],
1147  reversePointMap[e[1]]
1148  );
1149  }
1150  forAll(fem.edges(), i)
1151  {
1152  const edge& e = fem.edges()[i];
1153  newEdges[reverseFemEdgeMap[i]] = edge
1154  (
1155  reverseFemPointMap[e[0]],
1156  reverseFemPointMap[e[1]]
1157  );
1158  }
1159 
1160  pointField newEdgeDirections(newEdgeI);
1161  newEdgeDirections.rmap(edgeDirections(), reverseEdgeMap);
1162  newEdgeDirections.rmap(fem.edgeDirections(), reverseFemEdgeMap);
1163 
1164 
1165 
1166 
1167  // Normals
1168  // ~~~~~~~
1169 
1170  // Combine normals
1171  DynamicField<point> newNormals(normals().size()+fem.normals().size());
1172  newNormals.append(normals());
1173  newNormals.append(fem.normals());
1174 
1175 
1176  // Combine and re-index into newNormals
1177  labelListList newEdgeNormals(edgeNormals().size()+fem.edgeNormals().size());
1178  UIndirectList<labelList>(newEdgeNormals, reverseEdgeMap) =
1179  edgeNormals();
1180  UIndirectList<labelList>(newEdgeNormals, reverseFemEdgeMap) =
1181  fem.edgeNormals();
1182  forAll(reverseFemEdgeMap, i)
1183  {
1184  label mapI = reverseFemEdgeMap[i];
1185  labelList& en = newEdgeNormals[mapI];
1186  forAll(en, j)
1187  {
1188  en[j] += normals().size();
1189  }
1190  }
1191 
1192 
1193  // Combine and re-index into newFeaturePointNormals
1194  labelListList newFeaturePointNormals
1195  (
1196  featurePointNormals().size()
1197  + fem.featurePointNormals().size()
1198  );
1199 
1200  // Note: featurePointNormals only go up to nonFeatureStart
1202  (
1203  newFeaturePointNormals,
1204  SubList<label>(reversePointMap, featurePointNormals().size())
1205  ) = featurePointNormals();
1207  (
1208  newFeaturePointNormals,
1209  SubList<label>(reverseFemPointMap, fem.featurePointNormals().size())
1210  ) = fem.featurePointNormals();
1211  forAll(fem.featurePointNormals(), i)
1212  {
1213  label mapI = reverseFemPointMap[i];
1214  labelList& fn = newFeaturePointNormals[mapI];
1215  forAll(fn, j)
1216  {
1217  fn[j] += normals().size();
1218  }
1219  }
1220 
1221 
1222  // Combine regionEdges
1223  DynamicList<label> newRegionEdges
1224  (
1225  regionEdges().size()
1226  + fem.regionEdges().size()
1227  );
1228  forAll(regionEdges(), i)
1229  {
1230  newRegionEdges.append(reverseEdgeMap[regionEdges()[i]]);
1231  }
1232  forAll(fem.regionEdges(), i)
1233  {
1234  newRegionEdges.append(reverseFemEdgeMap[fem.regionEdges()[i]]);
1235  }
1236 
1237 
1238  // Assign
1239  // ~~~~~~
1240 
1241  // Transfer
1242  concaveStart_ = newConcaveStart;
1243  mixedStart_ = newMixedStart;
1244  nonFeatureStart_ = newNonFeatureStart;
1245 
1246  // Reset points and edges
1247  reset(move(newPoints), move(newEdges));
1248 
1249  // Transfer
1250  internalStart_ = newInternalStart;
1251  flatStart_ = newFlatStart;
1252  openStart_ = newOpenStart;
1253  multipleStart_ = newMultipleStart;
1254 
1255  edgeDirections_.transfer(newEdgeDirections);
1256 
1257  normals_.transfer(newNormals);
1258  edgeNormals_.transfer(newEdgeNormals);
1259  featurePointNormals_.transfer(newFeaturePointNormals);
1260 
1261  regionEdges_.transfer(newRegionEdges);
1262 
1263  pointTree_.clear();
1264  edgeTree_.clear();
1265  edgeTreesByType_.clear();
1266 }
1267 
1268 
1270 {
1271  // Points
1272  // ~~~~~~
1273 
1274  // From current points into new points
1275  labelList reversePointMap(identityMap(points().size()));
1276 
1277  // Flip convex and concave points
1278 
1279  label newPointi = 0;
1280  // Concave points become convex
1281  for (label i = concaveStart(); i < mixedStart(); i++)
1282  {
1283  reversePointMap[i] = newPointi++;
1284  }
1285  // Convex points become concave
1286  label newConcaveStart = newPointi;
1287  for (label i = 0; i < concaveStart(); i++)
1288  {
1289  reversePointMap[i] = newPointi++;
1290  }
1291 
1292 
1293  // Edges
1294  // ~~~~~~
1295 
1296  // From current edges into new edges
1297  labelList reverseEdgeMap(identityMap(edges().size()));
1298 
1299  // Flip external and internal edges
1300 
1301  label newEdgeI = 0;
1302  // Internal become external
1303  for (label i = internalStart(); i < flatStart(); i++)
1304  {
1305  reverseEdgeMap[i] = newEdgeI++;
1306  }
1307  // External become internal
1308  label newInternalStart = newEdgeI;
1309  for (label i = 0; i < internalStart(); i++)
1310  {
1311  reverseEdgeMap[i] = newEdgeI++;
1312  }
1313 
1314 
1315  pointField newPoints(points().size());
1316  newPoints.rmap(points(), reversePointMap);
1317 
1318  edgeList newEdges(edges().size());
1319  forAll(edges(), i)
1320  {
1321  const edge& e = edges()[i];
1322  newEdges[reverseEdgeMap[i]] = edge
1323  (
1324  reversePointMap[e[0]],
1325  reversePointMap[e[1]]
1326  );
1327  }
1328 
1329 
1330  // Normals are flipped
1331  // ~~~~~~~~~~~~~~~~~~~
1332 
1333  pointField newEdgeDirections(edges().size());
1334  newEdgeDirections.rmap(-1.0*edgeDirections(), reverseEdgeMap);
1335 
1336  pointField newNormals(-1.0*normals());
1337 
1338  labelListList newEdgeNormals(edgeNormals().size());
1339  UIndirectList<labelList>(newEdgeNormals, reverseEdgeMap) = edgeNormals();
1340 
1341  labelListList newFeaturePointNormals(featurePointNormals().size());
1342 
1343  // Note: featurePointNormals only go up to nonFeatureStart
1345  (
1346  newFeaturePointNormals,
1347  SubList<label>(reversePointMap, featurePointNormals().size())
1348  ) = featurePointNormals();
1349 
1350  labelList newRegionEdges(regionEdges().size());
1351  forAll(regionEdges(), i)
1352  {
1353  newRegionEdges[i] = reverseEdgeMap[regionEdges()[i]];
1354  }
1355 
1356  // Transfer
1357  concaveStart_ = newConcaveStart;
1358 
1359  // Reset points and edges
1360  reset(move(newPoints), move(newEdges));
1361 
1362  // Transfer
1363  internalStart_ = newInternalStart;
1364 
1365  edgeDirections_.transfer(newEdgeDirections);
1366  normals_.transfer(newNormals);
1367  edgeNormals_.transfer(newEdgeNormals);
1368  featurePointNormals_.transfer(newFeaturePointNormals);
1369  regionEdges_.transfer(newRegionEdges);
1370 
1371  pointTree_.clear();
1372  edgeTree_.clear();
1373  edgeTreesByType_.clear();
1374 }
1375 
1376 
1378 (
1379  const fileName& prefix,
1380  const bool verbose
1381 ) const
1382 {
1383  Info<< nl << "Writing extendedEdgeMesh components to " << prefix
1384  << endl;
1385 
1386  edgeMesh::write(prefix + "_edgeMesh.obj");
1387 
1388  if (!verbose) return;
1389 
1390  OBJstream convexFtPtStr(prefix + "_convexFeaturePts.obj");
1391  Info<< "Writing convex feature points to " << convexFtPtStr.name() << endl;
1392 
1393  for(label i = 0; i < concaveStart_; i++)
1394  {
1395  convexFtPtStr.write(points()[i]);
1396  }
1397 
1398  OBJstream concaveFtPtStr(prefix + "_concaveFeaturePts.obj");
1399  Info<< "Writing concave feature points to "
1400  << concaveFtPtStr.name() << endl;
1401 
1402  for(label i = concaveStart_; i < mixedStart_; i++)
1403  {
1404  convexFtPtStr.write(points()[i]);
1405  }
1406 
1407  OBJstream mixedFtPtStr(prefix + "_mixedFeaturePts.obj");
1408  Info<< "Writing mixed feature points to " << mixedFtPtStr.name() << endl;
1409 
1410  for(label i = mixedStart_; i < nonFeatureStart_; i++)
1411  {
1412  mixedFtPtStr.write(points()[i]);
1413  }
1414 
1415  OBJstream mixedFtPtStructureStr(prefix + "_mixedFeaturePtsStructure.obj");
1416  Info<< "Writing mixed feature point structure to "
1417  << mixedFtPtStructureStr.name() << endl;
1418 
1419  for(label i = mixedStart_; i < nonFeatureStart_; i++)
1420  {
1421  const labelList& ptEds = pointEdges()[i];
1422 
1423  forAll(ptEds, j)
1424  {
1425  const edge& e = edges()[ptEds[j]];
1426  mixedFtPtStructureStr.write
1427  (
1428  linePointRef(points()[e[0]],
1429  points()[e[1]])
1430  );
1431  }
1432  }
1433 
1434  OBJstream externalStr(prefix + "_externalEdges.obj");
1435  Info<< "Writing external edges to " << externalStr.name() << endl;
1436 
1437  for (label i = externalStart_; i < internalStart_; i++)
1438  {
1439  const edge& e = edges()[i];
1440  externalStr.write(linePointRef(points()[e[0]], points()[e[1]]));
1441  }
1442 
1443  OBJstream internalStr(prefix + "_internalEdges.obj");
1444  Info<< "Writing internal edges to " << internalStr.name() << endl;
1445 
1446  for (label i = internalStart_; i < flatStart_; i++)
1447  {
1448  const edge& e = edges()[i];
1449  internalStr.write(linePointRef(points()[e[0]], points()[e[1]]));
1450  }
1451 
1452  OBJstream flatStr(prefix + "_flatEdges.obj");
1453  Info<< "Writing flat edges to " << flatStr.name() << endl;
1454 
1455  for (label i = flatStart_; i < openStart_; i++)
1456  {
1457  const edge& e = edges()[i];
1458  flatStr.write(linePointRef(points()[e[0]], points()[e[1]]));
1459  }
1460 
1461  OBJstream openStr(prefix + "_openEdges.obj");
1462  Info<< "Writing open edges to " << openStr.name() << endl;
1463 
1464  for (label i = openStart_; i < multipleStart_; i++)
1465  {
1466  const edge& e = edges()[i];
1467  openStr.write(linePointRef(points()[e[0]], points()[e[1]]));
1468  }
1469 
1470  OBJstream multipleStr(prefix + "_multipleEdges.obj");
1471  Info<< "Writing multiple edges to " << multipleStr.name() << endl;
1472 
1473  for (label i = multipleStart_; i < edges().size(); i++)
1474  {
1475  const edge& e = edges()[i];
1476  multipleStr.write(linePointRef(points()[e[0]], points()[e[1]]));
1477  }
1478 
1479  OBJstream regionStr(prefix + "_regionEdges.obj");
1480  Info<< "Writing region edges to " << regionStr.name() << endl;
1481 
1482  forAll(regionEdges_, i)
1483  {
1484  const edge& e = edges()[regionEdges_[i]];
1485  regionStr.write(linePointRef(points()[e[0]], points()[e[1]]));
1486  }
1487 
1488  OBJstream edgeDirsStr(prefix + "_edgeDirections.obj");
1489  Info<< "Writing edge directions to " << edgeDirsStr.name() << endl;
1490 
1491  forAll(edgeDirections_, i)
1492  {
1493  const vector& eVec = edgeDirections_[i];
1494  const edge& e = edges()[i];
1495 
1496  edgeDirsStr.write
1497  (
1498  linePointRef(points()[e.start()], eVec + points()[e.start()])
1499  );
1500  }
1501 }
1502 
1503 
1505 {
1507 
1508  os << indent << "point classification :" << nl;
1509  os << incrIndent;
1510  os << indent << "convex feature points : "
1511  << setw(8) << concaveStart_-convexStart_
1512  << nl;
1513  os << indent << "concave feature points : "
1514  << setw(8) << mixedStart_-concaveStart_
1515  << nl;
1516  os << indent << "mixed feature points : "
1517  << setw(8) << nonFeatureStart_-mixedStart_
1518  << nl;
1519  os << indent << "other (non-feature) points : "
1520  << setw(8) << points().size()-nonFeatureStart_
1521  << nl;
1522  os << decrIndent;
1523 
1524  os << indent << "edge classification :" << nl;
1525  os << incrIndent;
1526  os << indent << "external (convex angle) edges : "
1527  << setw(8) << internalStart_-externalStart_
1528  << nl;
1529  os << indent << "internal (concave angle) edges : "
1530  << setw(8) << flatStart_-internalStart_
1531  << nl;
1532  os << indent << "flat region edges : "
1533  << setw(8) << openStart_-flatStart_
1534  << nl;
1535  os << indent << "open edges : "
1536  << setw(8) << multipleStart_-openStart_
1537  << nl;
1538  os << indent << "multiply connected edges : "
1539  << setw(8) << edges().size()-multipleStart_
1540  << nl;
1541  os << decrIndent;
1542 }
1543 
1544 
1545 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1546 
1547 Foam::Istream& Foam::operator>>
1548 (
1549  Istream& is,
1551 )
1552 {
1553  label type;
1554  is >> type;
1555 
1556  vt = static_cast<Foam::extendedEdgeMesh::sideVolumeType>(type);
1557 
1558  // Check state of Istream
1559  is.check("operator>>(Istream&, sideVolumeType&)");
1560 
1561  return is;
1562 }
1563 
1564 
1565 Foam::Ostream& Foam::operator<<
1566 (
1567  Ostream& os,
1569 )
1570 {
1571  os << static_cast<label>(vt);
1572 
1573  return os;
1574 }
1575 
1576 
1578 {
1579  // fileFormats::extendedEdgeMeshFormat::write(os, em.points_, em.edges_);
1580  os << "// points" << nl
1581  << em.points() << nl
1582  << "// edges" << nl
1583  << em.edges() << nl
1584  << "// concaveStart mixedStart nonFeatureStart" << nl
1585  << em.concaveStart_ << token::SPACE
1586  << em.mixedStart_ << token::SPACE
1587  << em.nonFeatureStart_ << nl
1588  << "// internalStart flatStart openStart multipleStart" << nl
1589  << em.internalStart_ << token::SPACE
1590  << em.flatStart_ << token::SPACE
1591  << em.openStart_ << token::SPACE
1592  << em.multipleStart_ << nl
1593  << "// normals" << nl
1594  << em.normals_ << nl
1595  << "// normal volume types" << nl
1596  << em.normalVolumeTypes_ << nl
1597  << "// normalDirections" << nl
1598  << em.normalDirections_ << nl
1599  << "// edgeNormals" << nl
1600  << em.edgeNormals_ << nl
1601  << "// featurePointNormals" << nl
1602  << em.featurePointNormals_ << nl
1603  << "// featurePointEdges" << nl
1604  << em.featurePointEdges_ << nl
1605  << "// regionEdges" << nl
1606  << em.regionEdges_
1607  << endl;
1608 
1609  // Check state of Ostream
1610  os.check("Ostream& operator<<(Ostream&, const extendedEdgeMesh&)");
1611 
1612  return os;
1613 }
1614 
1615 
1617 {
1618  // fileFormats::extendedEdgeMeshFormat::read(is, em.points_, em.edges_);
1619  is >> static_cast<edgeMesh&>(em)
1620  >> em.concaveStart_
1621  >> em.mixedStart_
1622  >> em.nonFeatureStart_
1623  >> em.internalStart_
1624  >> em.flatStart_
1625  >> em.openStart_
1626  >> em.multipleStart_
1627  >> em.normals_
1628  >> em.normalVolumeTypes_
1629  >> em.normalDirections_
1630  >> em.edgeNormals_
1631  >> em.featurePointNormals_
1632  >> em.featurePointEdges_
1633  >> em.regionEdges_;
1634 
1635  // Check state of Istream
1636  is.check("Istream& operator>>(Istream&, extendedEdgeMesh&)");
1637 
1638  return is;
1639 }
1640 
1641 
1642 // ************************************************************************* //
Istream and Ostream manipulators taking arguments.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Dynamically sized Field.
Definition: DynamicField.H:72
DynamicField< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:362
A HashTable with keys but without contents.
Definition: HashSet.H:62
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
OFstream which keeps track of vertices.
Definition: OBJstream.H:56
virtual Ostream & write(const char)
Write character.
Definition: OBJstream.C:81
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:120
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
bool hit() const
Is there a hit.
Definition: PointHit.H:120
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:54
const Point & hitPoint() const
Return hit point.
label index() const
Return index.
bool hit() const
Is there a hit.
A list of faces which address into the list of points.
const labelListList & edgeFaces() const
Return edge-face addressing.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A List obtained as a section of another List.
Definition: SubList.H:56
A List with indirect addressing.
Definition: UIndirectList.H:60
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
Points connected by edges.
Definition: edgeMesh.H:72
const edgeList & edges() const
Return edges.
Definition: edgeMeshI.H:68
void transfer(edgeMesh &)
Transfer the contents of the argument and annul the argument.
Definition: edgeMesh.C:196
virtual void writeStats(Ostream &) const
Definition: edgeMeshIO.C:117
const pointField & points() const
Return points.
Definition: edgeMeshI.H:62
virtual void clear()
Clear all storage.
Definition: edgeMesh.C:165
static void write(const fileName &, const edgeMesh &)
Write to file.
Definition: edgeMeshIO.C:87
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
bool connected(const edge &a) const
Return true if connected to given edge.
Definition: edgeI.H:103
Description of feature edges and points.
void allNearestFeatureEdges(const point &sample, const scalar searchRadiusSqr, pointIndexHitList &info) const
Find all the feature edges within searchDistSqr of sample.
void writeObj(const fileName &prefix, const bool verbose=true) const
Write all components of the extendedEdgeMesh as obj files.
labelList regionEdges_
Feature edges which are on the boundary between regions.
extendedEdgeMesh()
Construct null.
static const Foam::NamedEnum< sideVolumeType, 4 > sideVolumeTypeNames_
labelListList edgeNormals_
Indices of the normals that are adjacent to the feature edges.
static label nEdgeTypes
Number of possible feature edge types (i.e. number of slices)
void sortPointsAndEdges(const Patch &, const labelList &featureEdges, const labelList &regionFeatureEdges, const labelList &feaurePoints)
const indexedOctree< treeDataEdge > & edgeTree() const
Demand driven construction of octree for boundary edges.
label openStart() const
Return the index of the start of the open feature edges.
label nonFeatureStart() const
Return the index of the start of the non-feature points.
const vectorField & edgeDirections() const
Return the edgeDirection vectors.
static wordHashSet writeTypes()
pointStatus classifyFeaturePoint(label ptI) const
Classify the type of feature point. Requires valid stored member.
static const Foam::NamedEnum< edgeStatus, 6 > edgeStatusNames_
void nearestFeaturePoint(const point &sample, scalar searchDistSqr, pointIndexHit &info) const
Find nearest surface edge for the sample point.
label flatStart() const
Return the index of the start of the flat feature edges.
label nonFeatureStart_
Index of the start of the non-feature points.
const labelList & regionEdges() const
Return the feature edges which are on the boundary between.
const labelListList & edgeNormals() const
Return the indices of the normals that are adjacent to the.
vectorField edgeDirections_
Flat and open edges require the direction of the edge.
label openStart_
Index of the start of the open feature edges.
autoPtr< indexedOctree< treeDataEdge > > edgeTree_
Search tree for all edges.
labelListList featurePointEdges_
Indices of feature edges attached to feature points. The edges are.
label flatStart_
Index of the start of the flat feature edges.
static const Foam::NamedEnum< pointStatus, 4 > pointStatusNames_
static label nPointTypes
Number of possible point types (i.e. number of slices)
List< sideVolumeType > normalVolumeTypes_
Type per normal: which side of normal to mesh.
bool read(const fileName &, const word &ext)
Read from file. Chooses reader based on explicit extension.
static label externalStart_
Index of the start of the external feature edges - static as 0.
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.
scalar minDisconnectedDist(const pointIndexHitList &hitList) const
Return the minimum distance between disconnected edges.
static bool canWriteType(const word &ext, const bool verbose=false)
Can we write this file format type?
label internalStart_
Index of the start of the internal feature edges.
const vectorField & normals() const
Return the normals of the surfaces adjacent to the feature edges.
label multipleStart() const
Return the index of the start of the multiply-connected feature.
sideVolumeType
Normals point to the outside.
vectorField normals_
Normals of the features, to be referred to by index by both feature.
static label convexStart_
Index of the start of the convex feature points - static as 0.
label internalStart() const
Return the index of the start of the internal feature edges.
static bool canRead(const fileName &, const bool verbose=false)
Can we read this file format?
void allNearestFeaturePoints(const point &sample, scalar searchRadiusSqr, pointIndexHitList &info) const
Find all the feature points within searchDistSqr of sample.
static bool canReadType(const word &ext, const bool verbose=false)
Can we read this file format?
void nearestFeatureEdge(const point &sample, scalar searchDistSqr, pointIndexHit &info) const
Find nearest surface edge for the sample point.
virtual void writeStats(Ostream &os) const
Dump some information.
autoPtr< indexedOctree< treeDataPoint > > pointTree_
Search tree for all feature points.
virtual void clear()
Clear all storage.
label mixedStart_
Index of the start of the mixed type feature points.
void flipNormals()
Flip normals. All concave become convex, all internal external.
labelListList normalDirections_
Starting directions for the edges.
PtrList< indexedOctree< treeDataEdge > > edgeTreesByType_
Individual search trees for each type of edge.
~extendedEdgeMesh()
Destructor.
void transfer(extendedEdgeMesh &)
Transfer the contents of the argument and annul the argument.
const indexedOctree< treeDataPoint > & pointTree() const
Demand driven construction of octree for feature points.
const labelListList & featurePointNormals() const
Return the indices of the normals that are adjacent to the.
void add(const extendedEdgeMesh &)
Add extendedEdgeMesh. No filtering of duplicates.
label mixedStart() const
Return the index of the start of the mixed type feature points.
label multipleStart_
Index of the start of the multiply-connected feature edges.
void nearestFeatureEdgeByType(const point &sample, const scalarField &searchDistSqr, pointIndexHitList &info) const
Find the nearest point on each type of feature edge.
label concaveStart() const
Return the index of the start of the concave feature points.
labelListList featurePointNormals_
Indices of the normals that are adjacent to the feature points.
static scalar cosNormalAngleTol_
Angular closeness tolerance for treating normals as the same.
const PtrList< indexedOctree< treeDataEdge > > & edgeTreesByType() const
Demand driven construction of octree for boundary edges by type.
static wordHashSet readTypes()
label concaveStart_
Index of the start of the concave feature points.
A class for handling file names.
Definition: fileName.H:82
word ext() const
Return file name extension (part after last .)
Definition: fileName.C:299
Holds feature edges/points of surface.
const labelList & featureEdges() const
Return feature edge list.
const labelList & featurePoints() const
Return feature point list.
label nRegionEdges() const
Return number of region edges.
const triSurface & surface() const
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:54
Holds (reference to) pointField. Encapsulation of data needed for octree searches....
Definition: treeDataPoint.H:60
Triangulated surface description with patch information.
Definition: triSurface.H:69
A class for handling words, derived from string.
Definition: word.H:62
const pointField & points
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:105
scalar minDist(const List< pointIndexHit > &hitList)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
messageStream Info
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Istream & operator>>(Istream &, directionInfo &)
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:208
static const scalar GREAT
Definition: scalar.H:111
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Ostream & operator<<(Ostream &, const ensightPart &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
static const char nl
Definition: Ostream.H:260
dimensionedScalar cos(const dimensionedScalar &ds)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
label newPointi
Definition: readKivaGrid.H:495
scalarField samples(nIntervals, 0)