19#ifdef DEAL_II_WITH_OPENCASCADE
25# include <IGESControl_Controller.hxx>
26# include <IGESControl_Reader.hxx>
27# include <IGESControl_Writer.hxx>
28# include <STEPControl_Controller.hxx>
29# include <STEPControl_Reader.hxx>
30# include <STEPControl_Writer.hxx>
31# include <TopExp_Explorer.hxx>
33# include <TopoDS_Edge.hxx>
34# include <TopoDS_Face.hxx>
35# include <TopoDS_Shape.hxx>
39# if DEAL_II_OPENCASCADE_VERSION_GTE(7, 0, 0)
40# include <Standard_Transient.hxx>
42# include <Handle_Standard_Transient.hxx>
45# include <BRepAdaptor_Curve.hxx>
46# if DEAL_II_OPENCASCADE_VERSION_GTE(7, 6, 0)
47# include <BRepAlgoAPI_Section.hxx>
49# include <BRepAdaptor_HCompCurve.hxx>
50# include <BRepAdaptor_HCurve.hxx>
51# include <BRepAlgo_Section.hxx>
53# include <BRepAdaptor_Surface.hxx>
54# include <BRepBndLib.hxx>
55# include <BRepBuilderAPI_MakeEdge.hxx>
56# include <BRepBuilderAPI_Sewing.hxx>
57# include <BRepBuilderAPI_Transform.hxx>
58# include <BRepMesh_IncrementalMesh.hxx>
59# include <BRepTools.hxx>
60# include <BRep_Builder.hxx>
61# include <GCPnts_AbscissaPoint.hxx>
62# include <GeomAPI_Interpolate.hxx>
63# include <GeomAPI_ProjectPointOnCurve.hxx>
64# include <GeomAPI_ProjectPointOnSurf.hxx>
65# include <GeomConvert_CompCurveToBSplineCurve.hxx>
66# include <GeomLProp_SLProps.hxx>
67# include <Geom_BoundedCurve.hxx>
68# include <Geom_Plane.hxx>
69# include <IntCurvesFace_ShapeIntersector.hxx>
70# include <Poly_Triangulation.hxx>
71# include <ShapeAnalysis_Surface.hxx>
72# include <StlAPI_Reader.hxx>
73# include <StlAPI_Writer.hxx>
74# include <TColStd_HSequenceOfTransient.hxx>
75# include <TColStd_SequenceOfTransient.hxx>
76# include <TColgp_HArray1OfPnt.hxx>
77# include <TopLoc_Location.hxx>
90 std::tuple<unsigned int, unsigned int, unsigned int>
94 unsigned int n_faces = 0, n_edges = 0, n_vertices = 0;
95 for (exp.Init(shape, TopAbs_FACE); exp.More(); exp.Next(), ++n_faces)
98 for (exp.Init(shape, TopAbs_EDGE); exp.More(); exp.Next(), ++n_edges)
101 for (exp.Init(shape, TopAbs_VERTEX); exp.More(); exp.Next(), ++n_vertices)
104 return std::tuple<unsigned int, unsigned int, unsigned int>(n_faces,
111 std::vector<TopoDS_Face> &faces,
112 std::vector<TopoDS_Edge> &edges,
113 std::vector<TopoDS_Vertex> &
vertices)
120 for (exp.Init(shape, TopAbs_FACE); exp.More(); exp.Next())
122 faces.push_back(TopoDS::Face(exp.Current()));
124 for (exp.Init(shape, TopAbs_EDGE); exp.More(); exp.Next())
126 edges.push_back(TopoDS::Edge(exp.Current()));
128 for (exp.Init(shape, TopAbs_VERTEX); exp.More(); exp.Next())
130 vertices.push_back(TopoDS::Vertex(exp.Current()));
137 std::vector<TopoDS_Compound> &compounds,
138 std::vector<TopoDS_CompSolid> &compsolids,
139 std::vector<TopoDS_Solid> &solids,
140 std::vector<TopoDS_Shell> &shells,
141 std::vector<TopoDS_Wire> &wires)
144 compsolids.resize(0);
150 for (exp.Init(shape, TopAbs_COMPOUND); exp.More(); exp.Next())
152 compounds.push_back(TopoDS::Compound(exp.Current()));
154 for (exp.Init(shape, TopAbs_COMPSOLID); exp.More(); exp.Next())
156 compsolids.push_back(TopoDS::CompSolid(exp.Current()));
158 for (exp.Init(shape, TopAbs_SOLID); exp.More(); exp.Next())
160 solids.push_back(TopoDS::Solid(exp.Current()));
162 for (exp.Init(shape, TopAbs_SHELL); exp.More(); exp.Next())
164 shells.push_back(TopoDS::Shell(exp.Current()));
166 for (exp.Init(shape, TopAbs_WIRE); exp.More(); exp.Next())
168 wires.push_back(TopoDS::Wire(exp.Current()));
172 template <
int spacedim>
179 return gp_Pnt(p[0], 0, 0);
181 return gp_Pnt(p[0], p[1], 0);
183 return gp_Pnt(p[0], p[1], p[2]);
189 template <
int spacedim>
191 point(
const gp_Pnt &p,
const double tolerance)
199 "Cannot convert OpenCASCADE point to 1d if p.Y() != 0."));
202 "Cannot convert OpenCASCADE point to 1d if p.Z() != 0."));
207 "Cannot convert OpenCASCADE point to 2d if p.Z() != 0."));
221 const double tolerance)
223 const double rel_tol =
225 if (direction.
norm() > 0.0)
226 return (p1 * direction < p2 * direction - rel_tol);
228 for (
int d = dim; d >= 0; --d)
229 if (p1[d] < p2[d] - rel_tol)
231 else if (p2[d] < p1[d] - rel_tol)
241 read_IGES(
const std::string &filename,
const double scale_factor)
243 IGESControl_Reader reader;
244 IFSelect_ReturnStatus stat;
245 stat = reader.ReadFile(filename.c_str());
248 Standard_Boolean failsonly = Standard_False;
249 IFSelect_PrintCount mode = IFSelect_ItemsByEntity;
250 reader.PrintCheckLoad(failsonly, mode);
252 Standard_Integer nRoots = reader.TransferRoots();
261 scale.SetScale(Origin, scale_factor);
263 TopoDS_Shape sh = reader.OneShape();
264 BRepBuilderAPI_Transform trans(sh, scale);
266 return trans.Shape();
270 write_IGES(
const TopoDS_Shape &shape,
const std::string &filename)
272 IGESControl_Controller::Init();
273 IGESControl_Writer ICW(
"MM", 0);
274 Standard_Boolean ok = ICW.AddShape(shape);
277 Standard_Boolean OK = ICW.Write(filename.c_str());
285 StlAPI_Reader reader;
287 reader.Read(shape, filename.c_str());
294 const std::string &filename,
295 const double deflection,
296 const bool sew_different_faces,
297 const double sewer_tolerance,
298 const bool is_relative,
299 const double angular_deflection,
300 const bool in_parallel)
303 std::vector<TopoDS_Vertex>
vertices;
304 std::vector<TopoDS_Edge> edges;
305 std::vector<TopoDS_Face> faces;
307 const bool mesh_is_present =
308 std::none_of(faces.begin(), faces.end(), [&Loc](
const TopoDS_Face &face) {
309 Handle(Poly_Triangulation) theTriangulation =
310 BRep_Tool::Triangulation(face, Loc);
311 return theTriangulation.IsNull();
313 TopoDS_Shape shape_to_be_written = shape;
314 if (!mesh_is_present)
316 if (sew_different_faces)
318 BRepBuilderAPI_Sewing sewer(sewer_tolerance);
319 sewer.Add(shape_to_be_written);
321 shape_to_be_written = sewer.SewedShape();
324 shape_to_be_written = shape;
328 BRepMesh_IncrementalMesh mesh_im(shape_to_be_written,
335 StlAPI_Writer writer;
337# if DEAL_II_OPENCASCADE_VERSION_GTE(6, 9, 0)
339 const auto error = writer.Write(shape_to_be_written, filename.c_str());
342# if !DEAL_II_OPENCASCADE_VERSION_GTE(7, 2, 0)
353 writer.Write(shape_to_be_written, filename.c_str());
358 read_STEP(
const std::string &filename,
const double scale_factor)
360 STEPControl_Reader reader;
361 IFSelect_ReturnStatus stat;
362 stat = reader.ReadFile(filename.c_str());
365 Standard_Boolean failsonly = Standard_False;
366 IFSelect_PrintCount mode = IFSelect_ItemsByEntity;
367 reader.PrintCheckLoad(failsonly, mode);
369 Standard_Integer nRoots = reader.TransferRoots();
378 scale.SetScale(Origin, scale_factor);
380 TopoDS_Shape sh = reader.OneShape();
381 BRepBuilderAPI_Transform trans(sh, scale);
383 return trans.Shape();
387 write_STEP(
const TopoDS_Shape &shape,
const std::string &filename)
389 STEPControl_Controller::Init();
390 STEPControl_Writer SCW;
391 IFSelect_ReturnStatus status;
392 status = SCW.Transfer(shape, STEPControl_AsIs);
394 ExcMessage(
"Failed to add shape to STEP controller."));
396 status = SCW.Write(filename.c_str());
399 ExcMessage(
"Failed to write translated shape to STEP file."));
405 double tolerance = 0.0;
407 std::vector<TopoDS_Face> faces;
408 std::vector<TopoDS_Edge> edges;
409 std::vector<TopoDS_Vertex>
vertices;
414 tolerance = std::fmax(tolerance, BRep_Tool::Tolerance(vertex));
416 for (
const auto &edge : edges)
417 tolerance = std::fmax(tolerance, BRep_Tool::Tolerance(edge));
419 for (
const auto &face : faces)
420 tolerance = std::fmax(tolerance, BRep_Tool::Tolerance(face));
434 Handle(Geom_Plane) plane =
new Geom_Plane(c_x, c_y, c_z, c);
435# if DEAL_II_OPENCASCADE_VERSION_GTE(7, 6, 0)
436 BRepAlgoAPI_Section section(in_shape, plane);
438 BRepAlgo_Section section(in_shape, plane);
440 TopoDS_Shape edges = section.Shape();
445 join_edges(
const TopoDS_Shape &in_shape,
const double tolerance)
447 TopoDS_Edge out_shape;
448 const TopoDS_Shape &edges = in_shape;
449 std::vector<Handle_Geom_BoundedCurve> intersections;
453 gp_Pnt PIn(0.0, 0.0, 0.0);
454 gp_Pnt PFin(0.0, 0.0, 0.0);
455 gp_Pnt PMid(0.0, 0.0, 0.0);
456 TopExp_Explorer edgeExplorer(edges, TopAbs_EDGE);
458 while (edgeExplorer.More())
460 edge = TopoDS::Edge(edgeExplorer.Current());
461 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, L, First, Last);
462 intersections.push_back(Handle(Geom_BoundedCurve)::DownCast(curve));
468 unsigned int numIntersEdges = intersections.size();
471 GeomConvert_CompCurveToBSplineCurve convert_bspline(intersections[0]);
473 bool check =
false, one_added =
true, one_failed =
true;
474 std::vector<bool> added(numIntersEdges,
false);
476 while (one_added ==
true)
480 for (
unsigned int i = 1; i < numIntersEdges; ++i)
481 if (added[i] ==
false)
483 Handle(Geom_Curve) curve = intersections[i];
484 Handle(Geom_BoundedCurve) bcurve =
485 Handle(Geom_BoundedCurve)::DownCast(curve);
486 check = convert_bspline.Add(bcurve, tolerance,
false,
true, 0);
491 Handle(Geom_BoundedCurve) bcurve =
492 Handle(Geom_BoundedCurve)::DownCast(curve);
494 convert_bspline.Add(bcurve, tolerance,
false,
true, 0);
496 one_failed = one_failed || (check ==
false);
497 one_added = one_added || (check ==
true);
502 Assert(one_failed ==
false,
503 ExcMessage(
"Joining some of the Edges failed."));
505 Handle(Geom_Curve) bspline = convert_bspline.BSplineCurve();
507 out_shape = BRepBuilderAPI_MakeEdge(bspline);
516 const double tolerance)
520 gp_Pnt P0 =
point(origin);
523 dim > 1 ? direction[1] : 0,
524 dim > 2 ? direction[2] : 0));
528 gp_Pnt Pproj(0.0, 0.0, 0.0);
532 IntCurvesFace_ShapeIntersector Inters;
533 Inters.Load(in_shape, tolerance);
537 Inters.Perform(line, -RealLast(), +RealLast());
540 double minDistance = 1e7;
542 for (
int i = 0; i < Inters.NbPnt(); ++i)
544 const double distance =
point(origin).Distance(Inters.Pnt(i + 1));
547 if (distance < minDistance)
549 minDistance = distance;
562 const double tolerance)
564 unsigned int n_vertices = curve_points.size();
566 if (direction * direction > 0)
568 std::sort(curve_points.begin(),
571 return OpenCASCADE::point_compare(p1,
579 Handle(TColgp_HArray1OfPnt)
vertices =
580 new TColgp_HArray1OfPnt(1, n_vertices);
581 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
587 GeomAPI_Interpolate bspline_generator(
vertices, closed, tolerance);
588 bspline_generator.Perform();
589 Assert((bspline_generator.IsDone()),
590 ExcMessage(
"Interpolated bspline generation failed"));
592 Handle(Geom_BSplineCurve) bspline = bspline_generator.Curve();
593 TopoDS_Edge out_shape = BRepBuilderAPI_MakeEdge(bspline);
594 out_shape.Closed(closed);
600 template <
int spacedim>
601 std::vector<TopoDS_Edge>
609 std::map<unsigned int, std::pair<unsigned int, unsigned int>> vert_to_faces;
610 std::map<unsigned int, std::pair<unsigned int, unsigned int>> face_to_verts;
611 std::map<unsigned int, bool> visited_faces;
612 std::map<unsigned int, Point<spacedim>> vert_to_point;
614 unsigned int face_index;
618 if (cell->face(f)->at_boundary())
621 face_index = cell->face(f)->index();
622 const unsigned int v0 = cell->face(f)->vertex_index(0);
623 const unsigned int v1 = cell->face(f)->vertex_index(1);
624 face_to_verts[face_index].first =
v0;
625 face_to_verts[face_index].second =
v1;
626 visited_faces[face_index] =
false;
631 f, 0,
true,
false,
false)];
633 f, 1,
true,
false,
false)];
636 if (vert_to_faces.find(
v0) == vert_to_faces.end())
638 vert_to_faces[
v0].first = face_index;
642 vert_to_faces[
v0].second = face_index;
644 if (vert_to_faces.find(
v1) == vert_to_faces.end())
646 vert_to_faces[
v1].first = face_index;
650 vert_to_faces[
v1].second = face_index;
656 std::vector<TopoDS_Edge> interpolation_curves;
657 bool finished = (face_to_verts.empty());
658 face_index = finished ? 0 : face_to_verts.begin()->first;
660 while (finished ==
false)
662 const unsigned int start_point_index = face_to_verts[face_index].first;
663 unsigned int point_index = start_point_index;
666 std::vector<Point<spacedim>> pointlist;
669 visited_faces[face_index] =
true;
670 auto current_point = vert_to_point[point_index];
671 pointlist.push_back(current_point);
674 if (face_to_verts[face_index].
first != point_index)
675 point_index = face_to_verts[face_index].first;
677 point_index = face_to_verts[face_index].second;
680 if (vert_to_faces[point_index].
first != face_index)
681 face_index = vert_to_faces[point_index].first;
683 face_index = vert_to_faces[point_index].second;
685 while (point_index != start_point_index);
687 interpolation_curves.push_back(
691 for (
const auto &f : visited_faces)
692 if (f.second ==
false)
694 face_index = f.first;
699 return interpolation_curves;
704 std::tuple<Point<dim>, TopoDS_Shape, double,
double>
707 const double tolerance)
710 gp_Pnt Pproj =
point(origin);
712 double minDistance = 1e7;
713 gp_Pnt tmp_proj(0.0, 0.0, 0.0);
715 [[maybe_unused]]
unsigned int counter = 0;
716 unsigned int face_counter = 0;
718 TopoDS_Shape out_shape;
722 for (exp.Init(in_shape, TopAbs_FACE); exp.More(); exp.Next())
724 TopoDS_Face face = TopoDS::Face(exp.Current());
728 Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
730 ShapeAnalysis_Surface projector(SurfToProj);
731 gp_Pnt2d proj_params = projector.ValueOfUV(
point(origin), tolerance);
733 SurfToProj->D0(proj_params.X(), proj_params.Y(), tmp_proj);
736 if (distance < minDistance)
738 minDistance = distance;
752 if (face_counter == 0)
753 for (exp.Init(in_shape, TopAbs_EDGE); exp.More(); exp.Next())
755 TopoDS_Edge edge = TopoDS::Edge(exp.Current());
756 if (!BRep_Tool::Degenerated(edge))
764 Handle(Geom_Curve) CurveToProj =
765 BRep_Tool::Curve(edge, L, First, Last);
767 GeomAPI_ProjectPointOnCurve Proj(
point(origin), CurveToProj);
768 unsigned int num_proj_points = Proj.NbPoints();
769 if ((num_proj_points > 0) && (Proj.LowerDistance() < minDistance))
771 minDistance = Proj.LowerDistance();
772 Pproj = Proj.NearestPoint();
774 u = Proj.LowerDistanceParameter();
781 return std::tuple<Point<dim>, TopoDS_Shape, double,
double>(
790 const double tolerance)
792 std::tuple<Point<dim>, TopoDS_Shape, double,
double> ref =
794 return std::get<0>(ref);
800 const double tolerance)
803 std::tuple<Point<3>, TopoDS_Shape, double,
double> shape_and_params =
806 TopoDS_Shape &out_shape = std::get<1>(shape_and_params);
807 double &u = std::get<2>(shape_and_params);
808 double &v = std::get<3>(shape_and_params);
812 std::tuple<unsigned int, unsigned int, unsigned int>
numbers =
819 "Could not find normal: the shape containing the closest point has 0 faces."));
823 "Could not find normal: the shape containing the closest point has more than 1 face."));
827 exp.Init(out_shape, TopAbs_FACE);
828 TopoDS_Face face = TopoDS::Face(exp.Current());
834 push_forward(
const TopoDS_Shape &in_shape,
const double u,
const double v)
836 switch (in_shape.ShapeType())
840 BRepAdaptor_Surface surf(TopoDS::Face(in_shape));
845 BRepAdaptor_Curve curve(TopoDS::Edge(in_shape));
860 Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
861 GeomLProp_SLProps props(SurfToProj, u, v, 1, 1e-7);
862 gp_Pnt Value = props.Value();
864 gp_Dir Normal = props.Normal();
865 Assert(props.IsCurvatureDefined(),
866 ExcMessage(
"Curvature is not well defined!"));
867 Standard_Real Min_Curvature = props.MinCurvature();
868 Standard_Real Max_Curvature = props.MaxCurvature();
875 if (face.Orientation() == TopAbs_REVERSED)
890 template <
int spacedim>
895 BRepAdaptor_Surface surf(face);
896 const double u0 = surf.FirstUParameter();
897 const double u1 = surf.LastUParameter();
898 const double v0 = surf.FirstVParameter();
899 const double v1 = surf.LastVParameter();
901 std::vector<CellData<2>> cells;
902 std::vector<Point<spacedim>>
vertices;
911 for (
unsigned int i = 0; i < 4; ++i)
914 cells.push_back(cell);
918# include "utilities.inst"
Abstract base class for mapping classes.
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
numbers::NumberTraits< Number >::real_type norm() const
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcUnsupportedShape()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
TopoDS_Edge interpolation_curve(std::vector< Point< dim > > &curve_points, const Tensor< 1, dim > &direction=Tensor< 1, dim >(), const bool closed=false, const double tolerance=1e-7)
std::tuple< Point< 3 >, Tensor< 1, 3 >, double, double > push_forward_and_differential_forms(const TopoDS_Face &face, const double u, const double v, const double tolerance=1e-7)
void extract_compound_shapes(const TopoDS_Shape &shape, std::vector< TopoDS_Compound > &compounds, std::vector< TopoDS_CompSolid > &compsolids, std::vector< TopoDS_Solid > &solids, std::vector< TopoDS_Shell > &shells, std::vector< TopoDS_Wire > &wires)
Point< dim > push_forward(const TopoDS_Shape &in_shape, const double u, const double v)
std::tuple< unsigned int, unsigned int, unsigned int > count_elements(const TopoDS_Shape &shape)
void write_STEP(const TopoDS_Shape &shape, const std::string &filename)
std::tuple< Point< 3 >, Tensor< 1, 3 >, double, double > closest_point_and_differential_forms(const TopoDS_Shape &in_shape, const Point< 3 > &origin, const double tolerance=1e-7)
void extract_geometrical_shapes(const TopoDS_Shape &shape, std::vector< TopoDS_Face > &faces, std::vector< TopoDS_Edge > &edges, std::vector< TopoDS_Vertex > &vertices)
TopoDS_Edge join_edges(const TopoDS_Shape &in_shape, const double tolerance=1e-7)
TopoDS_Shape read_STEP(const std::string &filename, const double scale_factor=1e-3)
bool point_compare(const Point< dim > &p1, const Point< dim > &p2, const Tensor< 1, dim > &direction=Tensor< 1, dim >(), const double tolerance=1e-10)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Point< dim > closest_point(const TopoDS_Shape &in_shape, const Point< dim > &origin, const double tolerance=1e-7)
std::tuple< Point< dim >, TopoDS_Shape, double, double > project_point_and_pull_back(const TopoDS_Shape &in_shape, const Point< dim > &origin, const double tolerance=1e-7)
Point< dim > line_intersection(const TopoDS_Shape &in_shape, const Point< dim > &origin, const Tensor< 1, dim > &direction, const double tolerance=1e-7)
void write_STL(const TopoDS_Shape &shape, const std::string &filename, const double deflection, const bool sew_different_faces=false, const double sewer_tolerance=1e-6, const bool is_relative=false, const double angular_deflection=0.5, const bool in_parallel=false)
void create_triangulation(const TopoDS_Face &face, Triangulation< 2, spacedim > &tria)
TopoDS_Shape read_STL(const std::string &filename)
TopoDS_Shape intersect_plane(const TopoDS_Shape &in_shape, const double c_x, const double c_y, const double c_z, const double c, const double tolerance=1e-7)
std::vector< TopoDS_Edge > create_curves_from_triangulation_boundary(const Triangulation< 2, spacedim > &triangulation, const Mapping< 2, spacedim > &mapping=StaticMappingQ1< 2, spacedim >::mapping)
double get_shape_tolerance(const TopoDS_Shape &shape)
void write_IGES(const TopoDS_Shape &shape, const std::string &filename)
TopoDS_Shape read_IGES(const std::string &filename, const double scale_factor=1e-3)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::vector< unsigned int > vertices
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()