18 #ifdef DEAL_II_WITH_OPENCASCADE
24 # include <IGESControl_Controller.hxx>
25 # include <IGESControl_Reader.hxx>
26 # include <IGESControl_Writer.hxx>
27 # include <STEPControl_Controller.hxx>
28 # include <STEPControl_Reader.hxx>
29 # include <STEPControl_Writer.hxx>
30 # include <Standard_Version.hxx>
31 # include <TopExp_Explorer.hxx>
32 # include <TopoDS.hxx>
33 # include <TopoDS_Edge.hxx>
34 # include <TopoDS_Face.hxx>
35 # include <TopoDS_Shape.hxx>
40 # if (OCC_VERSION_MAJOR < 7)
41 # include <Handle_Standard_Transient.hxx>
43 # include <Standard_Transient.hxx>
46 # include <BRepAdaptor_Curve.hxx>
47 # include <BRepAdaptor_HCompCurve.hxx>
48 # include <BRepAdaptor_HCurve.hxx>
49 # include <BRepAdaptor_Surface.hxx>
50 # include <BRepAlgo_Section.hxx>
51 # include <BRepBndLib.hxx>
52 # include <BRepBuilderAPI_MakeEdge.hxx>
53 # include <BRepBuilderAPI_Sewing.hxx>
54 # include <BRepBuilderAPI_Transform.hxx>
55 # include <BRepMesh_IncrementalMesh.hxx>
56 # include <BRepTools.hxx>
57 # include <BRep_Builder.hxx>
58 # include <GCPnts_AbscissaPoint.hxx>
59 # include <GeomAPI_Interpolate.hxx>
60 # include <GeomAPI_ProjectPointOnCurve.hxx>
61 # include <GeomAPI_ProjectPointOnSurf.hxx>
62 # include <GeomConvert_CompCurveToBSplineCurve.hxx>
63 # include <GeomLProp_SLProps.hxx>
64 # include <Geom_BoundedCurve.hxx>
65 # include <Geom_Plane.hxx>
66 # include <IntCurvesFace_ShapeIntersector.hxx>
67 # include <Poly_Triangulation.hxx>
68 # include <ShapeAnalysis_Surface.hxx>
69 # include <StlAPI_Reader.hxx>
70 # include <StlAPI_Writer.hxx>
71 # include <TColStd_HSequenceOfTransient.hxx>
72 # include <TColStd_SequenceOfTransient.hxx>
73 # include <TColgp_HArray1OfPnt.hxx>
74 # include <TopLoc_Location.hxx>
75 # include <gp_Lin.hxx>
76 # include <gp_Pnt.hxx>
77 # include <gp_Vec.hxx>
87 std::tuple<unsigned int, unsigned int, unsigned int>
91 unsigned int n_faces = 0, n_edges = 0, n_vertices = 0;
92 for (exp.Init(shape, TopAbs_FACE); exp.More(); exp.Next(), ++n_faces)
95 for (exp.Init(shape, TopAbs_EDGE); exp.More(); exp.Next(), ++n_edges)
98 for (exp.Init(shape, TopAbs_VERTEX); exp.More(); exp.Next(), ++n_vertices)
101 return std::tuple<unsigned int, unsigned int, unsigned int>(n_faces,
108 std::vector<TopoDS_Face> & faces,
109 std::vector<TopoDS_Edge> & edges,
110 std::vector<TopoDS_Vertex> &
vertices)
117 for (exp.Init(shape, TopAbs_FACE); exp.More(); exp.Next())
119 faces.push_back(TopoDS::Face(exp.Current()));
121 for (exp.Init(shape, TopAbs_EDGE); exp.More(); exp.Next())
123 edges.push_back(TopoDS::Edge(exp.Current()));
125 for (exp.Init(shape, TopAbs_VERTEX); exp.More(); exp.Next())
134 std::vector<TopoDS_Compound> & compounds,
135 std::vector<TopoDS_CompSolid> &compsolids,
136 std::vector<TopoDS_Solid> & solids,
137 std::vector<TopoDS_Shell> & shells,
138 std::vector<TopoDS_Wire> & wires)
141 compsolids.resize(0);
147 for (exp.Init(shape, TopAbs_COMPOUND); exp.More(); exp.Next())
149 compounds.push_back(TopoDS::Compound(exp.Current()));
151 for (exp.Init(shape, TopAbs_COMPSOLID); exp.More(); exp.Next())
153 compsolids.push_back(TopoDS::CompSolid(exp.Current()));
155 for (exp.Init(shape, TopAbs_SOLID); exp.More(); exp.Next())
157 solids.push_back(TopoDS::Solid(exp.Current()));
159 for (exp.Init(shape, TopAbs_SHELL); exp.More(); exp.Next())
161 shells.push_back(TopoDS::Shell(exp.Current()));
163 for (exp.Init(shape, TopAbs_WIRE); exp.More(); exp.Next())
165 wires.push_back(TopoDS::Wire(exp.Current()));
169 template <
int spacedim>
176 return gp_Pnt(p[0], 0, 0);
178 return gp_Pnt(p[0], p[1], 0);
180 return gp_Pnt(p[0], p[1], p[2]);
186 template <
int spacedim>
188 point(
const gp_Pnt &p,
const double tolerance)
194 Assert(std::abs(p.Y()) <= tolerance,
196 "Cannot convert OpenCASCADE point to 1d if p.Y() != 0."));
197 Assert(std::abs(p.Z()) <= tolerance,
199 "Cannot convert OpenCASCADE point to 1d if p.Z() != 0."));
202 Assert(std::abs(p.Z()) <= tolerance,
204 "Cannot convert OpenCASCADE point to 2d if p.Z() != 0."));
218 const double tolerance)
220 const double rel_tol =
222 if (direction.
norm() > 0.0)
223 return (p1 * direction < p2 * direction - rel_tol);
225 for (
int d = dim;
d >= 0; --
d)
226 if (p1[
d] < p2[
d] - rel_tol)
228 else if (p2[
d] < p1[
d] - rel_tol)
238 read_IGES(
const std::string &filename,
const double scale_factor)
240 IGESControl_Reader reader;
241 IFSelect_ReturnStatus stat;
242 stat = reader.ReadFile(filename.c_str());
245 Standard_Boolean failsonly = Standard_False;
246 IFSelect_PrintCount mode = IFSelect_ItemsByEntity;
247 reader.PrintCheckLoad(failsonly, mode);
249 Standard_Integer nRoots = reader.TransferRoots();
258 scale.SetScale(Origin, scale_factor);
260 TopoDS_Shape sh = reader.OneShape();
261 BRepBuilderAPI_Transform trans(sh,
scale);
263 return trans.Shape();
267 write_IGES(
const TopoDS_Shape &shape,
const std::string &filename)
269 IGESControl_Controller::Init();
270 IGESControl_Writer ICW(
"MM", 0);
271 Standard_Boolean ok = ICW.AddShape(shape);
274 Standard_Boolean OK = ICW.Write(filename.c_str());
282 StlAPI_Reader reader;
284 reader.Read(shape, filename.c_str());
291 const std::string & filename,
292 const double deflection,
293 const bool sew_different_faces,
294 const double sewer_tolerance,
295 const bool is_relative,
296 const double angular_deflection,
297 const bool in_parallel)
300 std::vector<TopoDS_Vertex>
vertices;
301 std::vector<TopoDS_Edge> edges;
302 std::vector<TopoDS_Face> faces;
304 const bool mesh_is_present =
305 std::none_of(faces.begin(), faces.end(), [&Loc](
const TopoDS_Face &face) {
306 Handle(Poly_Triangulation) theTriangulation =
307 BRep_Tool::Triangulation(face, Loc);
308 return theTriangulation.IsNull();
310 TopoDS_Shape shape_to_be_written = shape;
311 if (!mesh_is_present)
313 if (sew_different_faces)
315 BRepBuilderAPI_Sewing sewer(sewer_tolerance);
316 sewer.Add(shape_to_be_written);
318 shape_to_be_written = sewer.SewedShape();
321 shape_to_be_written = shape;
325 BRepMesh_IncrementalMesh mesh_im(shape_to_be_written,
332 StlAPI_Writer writer;
334 # if ((OCC_VERSION_MAJOR * 100 + OCC_VERSION_MINOR * 10) >= 690)
336 const auto error = writer.Write(shape_to_be_written, filename.c_str());
339 # if ((OCC_VERSION_MAJOR * 100 + OCC_VERSION_MINOR * 10) < 720)
350 writer.Write(shape_to_be_written, filename.c_str());
355 read_STEP(
const std::string &filename,
const double scale_factor)
357 STEPControl_Reader reader;
358 IFSelect_ReturnStatus stat;
359 stat = reader.ReadFile(filename.c_str());
362 Standard_Boolean failsonly = Standard_False;
363 IFSelect_PrintCount mode = IFSelect_ItemsByEntity;
364 reader.PrintCheckLoad(failsonly, mode);
366 Standard_Integer nRoots = reader.TransferRoots();
375 scale.SetScale(Origin, scale_factor);
377 TopoDS_Shape sh = reader.OneShape();
378 BRepBuilderAPI_Transform trans(sh,
scale);
380 return trans.Shape();
384 write_STEP(
const TopoDS_Shape &shape,
const std::string &filename)
386 STEPControl_Controller::Init();
387 STEPControl_Writer SCW;
388 IFSelect_ReturnStatus status;
389 status = SCW.Transfer(shape, STEPControl_AsIs);
391 ExcMessage(
"Failed to add shape to STEP controller."));
393 status = SCW.Write(filename.c_str());
396 ExcMessage(
"Failed to write translated shape to STEP file."));
402 double tolerance = 0.0;
404 std::vector<TopoDS_Face> faces;
405 std::vector<TopoDS_Edge> edges;
406 std::vector<TopoDS_Vertex>
vertices;
411 tolerance = std::fmax(tolerance, BRep_Tool::Tolerance(vertex));
413 for (
const auto &edge : edges)
414 tolerance = std::fmax(tolerance, BRep_Tool::Tolerance(edge));
416 for (
const auto &face : faces)
417 tolerance = std::fmax(tolerance, BRep_Tool::Tolerance(face));
431 Handle(Geom_Plane) plane =
new Geom_Plane(c_x, c_y, c_z, c);
432 BRepAlgo_Section section(in_shape, plane);
433 TopoDS_Shape edges = section.Shape();
438 join_edges(
const TopoDS_Shape &in_shape,
const double tolerance)
440 TopoDS_Edge out_shape;
441 const TopoDS_Shape & edges = in_shape;
442 std::vector<Handle_Geom_BoundedCurve> intersections;
446 gp_Pnt PIn(0.0, 0.0, 0.0);
447 gp_Pnt PFin(0.0, 0.0, 0.0);
448 gp_Pnt PMid(0.0, 0.0, 0.0);
449 TopExp_Explorer edgeExplorer(edges, TopAbs_EDGE);
451 while (edgeExplorer.More())
453 edge = TopoDS::Edge(edgeExplorer.Current());
454 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge,
L, First, Last);
455 intersections.push_back(Handle(Geom_BoundedCurve)::DownCast(curve));
461 unsigned int numIntersEdges = intersections.size();
464 GeomConvert_CompCurveToBSplineCurve convert_bspline(intersections[0]);
466 bool check =
false, one_added =
true, one_failed =
true;
467 std::vector<bool> added(numIntersEdges,
false);
469 while (one_added ==
true)
473 for (
unsigned int i = 1; i < numIntersEdges; ++i)
474 if (added[i] ==
false)
476 Handle(Geom_Curve) curve = intersections[i];
477 Handle(Geom_BoundedCurve) bcurve =
478 Handle(Geom_BoundedCurve)::DownCast(curve);
479 check = convert_bspline.Add(bcurve, tolerance,
false,
true, 0);
484 Handle(Geom_BoundedCurve) bcurve =
485 Handle(Geom_BoundedCurve)::DownCast(curve);
487 convert_bspline.Add(bcurve, tolerance,
false,
true, 0);
489 one_failed = one_failed || (check ==
false);
490 one_added = one_added || (check ==
true);
495 Assert(one_failed ==
false,
496 ExcMessage(
"Joining some of the Edges failed."));
498 Handle(Geom_Curve) bspline = convert_bspline.BSplineCurve();
500 out_shape = BRepBuilderAPI_MakeEdge(bspline);
509 const double tolerance)
513 gp_Pnt P0 =
point(origin);
516 dim > 1 ? direction[1] : 0,
517 dim > 2 ? direction[2] : 0));
521 gp_Pnt Pproj(0.0, 0.0, 0.0);
525 IntCurvesFace_ShapeIntersector Inters;
526 Inters.Load(in_shape, tolerance);
530 Inters.Perform(line, -RealLast(), +RealLast());
533 double minDistance = 1e7;
535 for (
int i = 0; i < Inters.NbPnt(); ++i)
537 const double distance =
point(origin).Distance(Inters.Pnt(i + 1));
540 if (distance < minDistance)
542 minDistance = distance;
543 result = point<dim>(Inters.Pnt(i + 1));
555 const double tolerance)
557 unsigned int n_vertices = curve_points.size();
559 if (direction * direction > 0)
561 std::sort(curve_points.begin(),
564 return OpenCASCADE::point_compare(p1,
572 Handle(TColgp_HArray1OfPnt)
vertices =
573 new TColgp_HArray1OfPnt(1, n_vertices);
574 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
580 GeomAPI_Interpolate bspline_generator(
vertices, closed, tolerance);
581 bspline_generator.Perform();
582 Assert((bspline_generator.IsDone()),
583 ExcMessage(
"Interpolated bspline generation failed"));
585 Handle(Geom_BSplineCurve) bspline = bspline_generator.Curve();
586 TopoDS_Edge out_shape = BRepBuilderAPI_MakeEdge(bspline);
587 out_shape.Closed(closed);
593 template <
int spacedim>
594 std::vector<TopoDS_Edge>
602 std::map<unsigned int, std::pair<unsigned int, unsigned int>> vert_to_faces;
603 std::map<unsigned int, std::pair<unsigned int, unsigned int>> face_to_verts;
604 std::map<unsigned int, bool> visited_faces;
605 std::map<unsigned int, Point<spacedim>> vert_to_point;
607 unsigned int face_index;
609 for (
const auto &cell :
triangulation.active_cell_iterators())
611 if (cell->face(f)->at_boundary())
614 face_index = cell->face(f)->index();
615 const unsigned int v0 = cell->face(f)->vertex_index(0);
616 const unsigned int v1 = cell->face(f)->vertex_index(1);
617 face_to_verts[face_index].first = v0;
618 face_to_verts[face_index].second = v1;
619 visited_faces[face_index] =
false;
625 f, 0,
true,
false,
false)];
627 f, 1,
true,
false,
false)];
630 if (vert_to_faces.find(v0) == vert_to_faces.end())
632 vert_to_faces[v0].first = face_index;
636 vert_to_faces[v0].second = face_index;
638 if (vert_to_faces.find(v1) == vert_to_faces.end())
640 vert_to_faces[v1].first = face_index;
644 vert_to_faces[v1].second = face_index;
650 std::vector<TopoDS_Edge> interpolation_curves;
651 bool finished = (face_to_verts.size() == 0);
652 face_index = finished ? 0 : face_to_verts.begin()->first;
654 while (finished ==
false)
656 const unsigned int start_point_index = face_to_verts[face_index].first;
657 unsigned int point_index = start_point_index;
660 std::vector<Point<spacedim>> pointlist;
663 visited_faces[face_index] =
true;
664 auto current_point = vert_to_point[point_index];
665 pointlist.push_back(current_point);
668 if (face_to_verts[face_index].
first != point_index)
669 point_index = face_to_verts[face_index].first;
671 point_index = face_to_verts[face_index].second;
674 if (vert_to_faces[point_index].
first != face_index)
675 face_index = vert_to_faces[point_index].first;
677 face_index = vert_to_faces[point_index].second;
679 while (point_index != start_point_index);
681 interpolation_curves.push_back(
685 for (
const auto &f : visited_faces)
686 if (f.second ==
false)
688 face_index = f.first;
693 return interpolation_curves;
698 std::tuple<Point<dim>, TopoDS_Shape,
double,
double>
701 const double tolerance)
704 gp_Pnt Pproj =
point(origin);
706 double minDistance = 1e7;
707 gp_Pnt tmp_proj(0.0, 0.0, 0.0);
709 unsigned int counter = 0;
710 unsigned int face_counter = 0;
712 TopoDS_Shape out_shape;
716 for (exp.Init(in_shape, TopAbs_FACE); exp.More(); exp.Next())
718 TopoDS_Face face = TopoDS::Face(exp.Current());
722 Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
724 ShapeAnalysis_Surface projector(SurfToProj);
725 gp_Pnt2d proj_params = projector.ValueOfUV(
point(origin), tolerance);
727 SurfToProj->D0(proj_params.X(), proj_params.Y(), tmp_proj);
729 double distance = point<dim>(tmp_proj).distance(origin);
730 if (distance < minDistance)
732 minDistance = distance;
746 if (face_counter == 0)
747 for (exp.Init(in_shape, TopAbs_EDGE); exp.More(); exp.Next())
749 TopoDS_Edge edge = TopoDS::Edge(exp.Current());
750 if (!BRep_Tool::Degenerated(edge))
758 Handle(Geom_Curve) CurveToProj =
759 BRep_Tool::Curve(edge,
L, First, Last);
761 GeomAPI_ProjectPointOnCurve Proj(
point(origin), CurveToProj);
762 unsigned int num_proj_points = Proj.NbPoints();
763 if ((num_proj_points > 0) && (Proj.LowerDistance() < minDistance))
765 minDistance = Proj.LowerDistance();
766 Pproj = Proj.NearestPoint();
768 u = Proj.LowerDistanceParameter();
775 return std::tuple<Point<dim>, TopoDS_Shape,
double,
double>(
776 point<dim>(Pproj), out_shape, u, v);
784 const double tolerance)
786 std::tuple<Point<dim>, TopoDS_Shape,
double,
double> ref =
788 return std::get<0>(ref);
794 const double tolerance)
797 std::tuple<Point<3>, TopoDS_Shape,
double,
double> shape_and_params =
800 TopoDS_Shape &out_shape = std::get<1>(shape_and_params);
801 double & u = std::get<2>(shape_and_params);
802 double & v = std::get<3>(shape_and_params);
806 std::tuple<unsigned int, unsigned int, unsigned int>
numbers =
813 "Could not find normal: the shape containing the closest point has 0 faces."));
817 "Could not find normal: the shape containing the closest point has more than 1 face."));
821 exp.Init(out_shape, TopAbs_FACE);
822 TopoDS_Face face = TopoDS::Face(exp.Current());
828 push_forward(
const TopoDS_Shape &in_shape,
const double u,
const double v)
830 switch (in_shape.ShapeType())
834 BRepAdaptor_Surface surf(TopoDS::Face(in_shape));
835 return point<dim>(surf.Value(u, v));
839 BRepAdaptor_Curve curve(TopoDS::Edge(in_shape));
840 return point<dim>(curve.Value(u));
854 Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
855 GeomLProp_SLProps props(SurfToProj, u, v, 1, 1
e-7);
856 gp_Pnt Value = props.Value();
858 gp_Dir Normal = props.Normal();
859 Assert(props.IsCurvatureDefined(),
860 ExcMessage(
"Curvature is not well defined!"));
861 Standard_Real Min_Curvature = props.MinCurvature();
862 Standard_Real Max_Curvature = props.MaxCurvature();
869 if (face.Orientation() == TopAbs_REVERSED)
884 template <
int spacedim>
889 BRepAdaptor_Surface surf(face);
890 const double u0 = surf.FirstUParameter();
891 const double u1 = surf.LastUParameter();
892 const double v0 = surf.FirstVParameter();
893 const double v1 = surf.LastVParameter();
895 std::vector<CellData<2>> cells;
896 std::vector<Point<spacedim>>
vertices;
899 vertices.push_back(point<spacedim>(surf.Value(u0, v0)));
900 vertices.push_back(point<spacedim>(surf.Value(u1, v0)));
901 vertices.push_back(point<spacedim>(surf.Value(u0, v1)));
902 vertices.push_back(point<spacedim>(surf.Value(u1, v1)));
905 for (
unsigned int i = 0; i < 4; ++i)
908 cells.push_back(cell);
912 # include "utilities.inst"