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