16 #include <deal.II/opencascade/utilities.h> 18 #ifdef DEAL_II_WITH_OPENCASCADE 20 #include <deal.II/base/point.h> 21 #include <deal.II/base/utilities.h> 22 #include <deal.II/base/exceptions.h> 24 #include <boost/bind.hpp> 30 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
32 #include <IGESControl_Controller.hxx> 33 #include <IGESControl_Reader.hxx> 34 #include <IGESControl_Writer.hxx> 36 #include <STEPControl_Controller.hxx> 37 #include <STEPControl_Reader.hxx> 38 #include <STEPControl_Writer.hxx> 41 #include <TopoDS_Shape.hxx> 42 #include <TopoDS_Face.hxx> 43 #include <TopoDS_Edge.hxx> 44 #include <TopExp_Explorer.hxx> 46 #include <Standard_Version.hxx> 47 #if (OCC_VERSION_MAJOR < 7) 48 # include <Handle_Standard_Transient.hxx> 50 # include <Standard_Transient.hxx> 53 #include <TColStd_SequenceOfTransient.hxx> 54 #include <TColStd_HSequenceOfTransient.hxx> 55 #include <TColgp_HArray1OfPnt.hxx> 60 #include <GeomAPI_ProjectPointOnSurf.hxx> 61 #include <GeomAPI_ProjectPointOnCurve.hxx> 62 #include <IntCurvesFace_ShapeIntersector.hxx> 64 #include <BRepTools.hxx> 65 #include <BRepAdaptor_Curve.hxx> 66 #include <BRepAdaptor_HCurve.hxx> 67 #include <BRepAdaptor_HCompCurve.hxx> 68 #include <BRepAdaptor_Surface.hxx> 69 #include <BRep_Builder.hxx> 70 #include <BRepBuilderAPI_Transform.hxx> 71 #include <BRepBuilderAPI_MakeEdge.hxx> 72 #include <BRepAlgo_Section.hxx> 74 #include <Geom_Plane.hxx> 75 #include <Geom_BoundedCurve.hxx> 76 #include <GeomAPI_Interpolate.hxx> 77 #include <GeomConvert_CompCurveToBSplineCurve.hxx> 78 #include <GeomLProp_SLProps.hxx> 80 #include <GCPnts_AbscissaPoint.hxx> 81 #include <ShapeAnalysis_Surface.hxx> 83 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
88 DEAL_II_NAMESPACE_OPEN
92 std_cxx11::tuple<unsigned int, unsigned int, unsigned int>
96 unsigned int n_faces=0, n_edges=0, n_vertices=0;
97 for (exp.Init(shape, TopAbs_FACE);
98 exp.More(); exp.Next(), ++n_faces)
100 for (exp.Init(shape, TopAbs_EDGE);
101 exp.More(); exp.Next(), ++n_edges)
103 for (exp.Init(shape, TopAbs_VERTEX);
104 exp.More(); exp.Next(), ++n_vertices)
106 return std_cxx11::tuple<unsigned int, unsigned int, unsigned int>(n_faces, n_edges, n_vertices);
110 std::vector<TopoDS_Face> &faces,
111 std::vector<TopoDS_Edge> &edges,
112 std::vector<TopoDS_Vertex> &vertices)
119 for (exp.Init(shape, TopAbs_FACE); exp.More(); exp.Next())
121 faces.push_back(TopoDS::Face(exp.Current()));
123 for (exp.Init(shape, TopAbs_EDGE); exp.More(); exp.Next())
125 edges.push_back(TopoDS::Edge(exp.Current()));
127 for (exp.Init(shape, TopAbs_VERTEX); exp.More(); exp.Next())
129 vertices.push_back(TopoDS::Vertex(exp.Current()));
135 std::vector<TopoDS_Compound> &compounds,
136 std::vector<TopoDS_CompSolid> &compsolids,
137 std::vector<TopoDS_Solid> &solids,
138 std::vector<TopoDS_Shell> &shells,
139 std::vector<TopoDS_Wire> &wires)
142 compsolids.resize(0);
148 for (exp.Init(shape, TopAbs_COMPOUND); exp.More(); exp.Next())
150 compounds.push_back(TopoDS::Compound(exp.Current()));
152 for (exp.Init(shape, TopAbs_COMPSOLID); exp.More(); exp.Next())
154 compsolids.push_back(TopoDS::CompSolid(exp.Current()));
156 for (exp.Init(shape, TopAbs_SOLID); exp.More(); exp.Next())
158 solids.push_back(TopoDS::Solid(exp.Current()));
160 for (exp.Init(shape, TopAbs_SHELL); exp.More(); exp.Next())
162 shells.push_back(TopoDS::Shell(exp.Current()));
164 for (exp.Init(shape, TopAbs_WIRE); exp.More(); exp.Next())
166 wires.push_back(TopoDS::Wire(exp.Current()));
172 return gp_Pnt(p(0), p(1), p(2));
178 return Point<3>(p.X(), p.Y(), p.Z());
183 const double tolerance)
185 const double rel_tol=std::max(tolerance, std::max(p1.
norm(), p2.
norm())*tolerance);
186 if (direction.
norm() > 0.0)
187 return (p1*direction < p2*direction-rel_tol);
189 for (
int d=2; d>=0; --d)
190 if (p1[d] < p2[d]-rel_tol)
192 else if (p2[d] < p1[d]-rel_tol)
202 const double scale_factor)
204 IGESControl_Reader reader;
205 IFSelect_ReturnStatus stat;
206 stat = reader.ReadFile(filename.c_str());
210 Standard_Boolean failsonly = Standard_False;
211 IFSelect_PrintCount mode = IFSelect_ItemsByEntity;
212 reader.PrintCheckLoad (failsonly, mode);
214 Standard_Integer nRoots = reader.TransferRoots();
224 scale.SetScale (Origin, scale_factor);
226 TopoDS_Shape sh = reader.OneShape();
227 BRepBuilderAPI_Transform trans(sh, scale);
229 return trans.Shape();
233 const std::string &filename)
235 IGESControl_Controller::Init();
236 IGESControl_Writer ICW (
"MM", 0);
237 Standard_Boolean ok = ICW.AddShape (shape);
240 Standard_Boolean OK = ICW.Write (filename.c_str());
245 const double scale_factor)
247 STEPControl_Reader reader;
248 IFSelect_ReturnStatus stat;
249 stat = reader.ReadFile(filename.c_str());
253 Standard_Boolean failsonly = Standard_False;
254 IFSelect_PrintCount mode = IFSelect_ItemsByEntity;
255 reader.PrintCheckLoad (failsonly, mode);
257 Standard_Integer nRoots = reader.TransferRoots();
267 scale.SetScale (Origin, scale_factor);
269 TopoDS_Shape sh = reader.OneShape();
270 BRepBuilderAPI_Transform trans(sh, scale);
272 return trans.Shape();
276 const std::string &filename)
278 STEPControl_Controller::Init();
279 STEPControl_Writer SCW;
280 IFSelect_ReturnStatus status;
281 status = SCW.Transfer(shape, STEPControl_AsIs);
284 status = SCW.Write(filename.c_str());
286 AssertThrow(status == IFSelect_RetDone,
ExcMessage(
"Failed to write translated shape to STEP file."));
291 double tolerance = 0.0;
293 std::vector<TopoDS_Face> faces;
294 std::vector<TopoDS_Edge> edges;
295 std::vector<TopoDS_Vertex> vertices;
302 for (
unsigned int i=0; i<vertices.size(); ++i)
303 tolerance = fmax(tolerance,BRep_Tool::Tolerance(vertices[i]));
305 for (
unsigned int i=0; i<edges.size(); ++i)
306 tolerance = fmax(tolerance,BRep_Tool::Tolerance(edges[i]));
308 for (
unsigned int i=0; i<faces.size(); ++i)
309 tolerance = fmax(tolerance,BRep_Tool::Tolerance(faces[i]));
322 Handle(Geom_Plane) plane =
new Geom_Plane(c_x,c_y,c_z,c);
323 BRepAlgo_Section section(in_shape, plane);
324 TopoDS_Shape edges = section.Shape();
329 const double tolerance)
331 TopoDS_Edge out_shape;
332 TopoDS_Shape edges = in_shape;
333 std::vector<Handle_Geom_BoundedCurve> intersections;
337 gp_Pnt PIn(0.0,0.0,0.0);
338 gp_Pnt PFin(0.0,0.0,0.0);
339 gp_Pnt PMid(0.0,0.0,0.0);
340 TopExp_Explorer edgeExplorer(edges , TopAbs_EDGE);
342 while (edgeExplorer.More())
344 edge = TopoDS::Edge(edgeExplorer.Current());
345 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge,L,First,Last);
346 intersections.push_back(Handle(Geom_BoundedCurve)::DownCast(curve));
352 unsigned int numIntersEdges = intersections.size();
355 GeomConvert_CompCurveToBSplineCurve convert_bspline(intersections[0]);
357 bool check =
false, one_added =
true, one_failed=
true;
358 std::vector<bool> added(numIntersEdges,
false);
360 while (one_added ==
true)
364 for (
unsigned int i=1; i<numIntersEdges; ++i)
365 if (added[i] ==
false)
367 Handle(Geom_Curve) curve = intersections[i];
368 Handle(Geom_BoundedCurve) bcurve = Handle(Geom_BoundedCurve)::DownCast(curve);
369 check = convert_bspline.Add(bcurve,tolerance,0,1,0);
373 Handle(Geom_BoundedCurve) bcurve = Handle(Geom_BoundedCurve)::DownCast(curve);
374 check = convert_bspline.Add(bcurve,tolerance,0,1,0);
376 one_failed = one_failed || (check ==
false);
377 one_added = one_added || (check ==
true);
382 Assert(one_failed ==
false,
383 ExcMessage(
"Joining some of the Edges failed."));
385 Handle(Geom_Curve) bspline = convert_bspline.BSplineCurve();
387 out_shape = BRepBuilderAPI_MakeEdge(bspline);
395 const double tolerance)
399 gp_Pnt P0 =
point(origin);
400 gp_Ax1 gpaxis(P0, gp_Dir(direction[0], direction[1], direction[2]));
404 gp_Pnt Pproj(0.0,0.0,0.0);
408 IntCurvesFace_ShapeIntersector Inters;
409 Inters.Load(in_shape,tolerance);
413 Inters.Perform(line,-RealLast(),+RealLast());
416 double minDistance = 1e7;
418 for (
int i=0; i<Inters.NbPnt(); ++i)
420 const double distance =
point(origin).Distance(Inters.Pnt(i+1));
422 if (distance < minDistance)
424 minDistance = distance;
425 result =
point(Inters.Pnt(i+1));
435 const double tolerance)
438 unsigned int n_vertices = curve_points.size();
440 if (direction*direction > 0)
442 std::sort(curve_points.begin(), curve_points.end(),
447 Handle(TColgp_HArray1OfPnt) vertices =
new TColgp_HArray1OfPnt(1,n_vertices);
448 for (
unsigned int vertex=0; vertex<n_vertices; ++vertex)
450 vertices->SetValue(vertex+1,
point(curve_points[vertex]));
454 GeomAPI_Interpolate bspline_generator(vertices, closed, tolerance);
455 bspline_generator.Perform();
456 Assert( (bspline_generator.IsDone()),
ExcMessage(
"Interpolated bspline generation failed"));
458 Handle(Geom_BSplineCurve) bspline = bspline_generator.Curve();
459 TopoDS_Edge out_shape = BRepBuilderAPI_MakeEdge(bspline);
463 std_cxx11::tuple<Point<3>, TopoDS_Shape, double,
double>
466 const double tolerance)
469 gp_Pnt Pproj =
point(origin);
471 double minDistance = 1e7;
472 gp_Pnt tmp_proj(0.0,0.0,0.0);
474 unsigned int counter = 0;
475 unsigned int face_counter = 0;
477 TopoDS_Shape out_shape;
481 for (exp.Init(in_shape, TopAbs_FACE); exp.More(); exp.Next())
483 TopoDS_Face face = TopoDS::Face(exp.Current());
487 Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
489 ShapeAnalysis_Surface projector(SurfToProj);
490 gp_Pnt2d proj_params = projector.ValueOfUV(
point(origin), tolerance);
492 SurfToProj->D0(proj_params.X(),proj_params.Y(),tmp_proj);
495 if (distance < minDistance)
497 minDistance = distance;
512 for (exp.Init(in_shape, TopAbs_EDGE); exp.More(); exp.Next())
514 TopoDS_Edge edge = TopoDS::Edge(exp.Current());
515 if (!BRep_Tool::Degenerated(edge))
523 Handle(Geom_Curve) CurveToProj = BRep_Tool::Curve(edge,L,First,Last);
525 GeomAPI_ProjectPointOnCurve Proj(
point(origin),CurveToProj);
526 unsigned int num_proj_points = Proj.NbPoints();
527 if ((num_proj_points > 0) && (Proj.LowerDistance() < minDistance))
529 minDistance = Proj.LowerDistance();
530 Pproj = Proj.NearestPoint();
532 u=Proj.LowerDistanceParameter();
539 return std_cxx11::tuple<Point<3>, TopoDS_Shape, double,
double>
540 (
point(Pproj),out_shape, u, v);
546 const double tolerance)
548 std_cxx11::tuple<Point<3>, TopoDS_Shape, double,
double>
550 return std_cxx11::get<0>(ref);
553 std_cxx11::tuple<Point<3>,
Tensor<1,3>, double,
double>
556 const double tolerance)
559 std_cxx11::tuple<Point<3>, TopoDS_Shape, double,
double>
564 TopoDS_Shape &out_shape = std_cxx11::get<1>(shape_and_params);
565 double &u = std_cxx11::get<2>(shape_and_params);
566 double &v = std_cxx11::get<3>(shape_and_params);
570 std_cxx11::tuple<unsigned int, unsigned int, unsigned int>
numbers =
575 ExcMessage(
"Could not find normal: the shape containing the closest point has 0 faces."));
577 ExcMessage(
"Could not find normal: the shape containing the closest point has more than 1 face."));
581 exp.Init(out_shape, TopAbs_FACE);
582 TopoDS_Face face = TopoDS::Face(exp.Current());
590 switch (in_shape.ShapeType())
594 BRepAdaptor_Surface surf(TopoDS::Face(in_shape));
595 return point(surf.Value(u,v));
599 BRepAdaptor_Curve curve(TopoDS::Edge(in_shape));
600 return point(curve.Value(u));
608 std_cxx11::tuple<Point<3>,
Tensor<1,3>, double,
double>
614 Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
615 GeomLProp_SLProps props(SurfToProj, u, v, 1, 1e-7);
616 gp_Pnt Value = props.Value();
618 gp_Dir Normal = props.Normal();
619 Assert(props.IsCurvatureDefined(),
ExcMessage(
"Curvature is not well defined!"));
620 Standard_Real Min_Curvature = props.MinCurvature();
621 Standard_Real Max_Curvature = props.MaxCurvature();
628 if (face.Orientation()==TopAbs_REVERSED)
635 return std_cxx11::tuple<Point<3>,
Tensor<1,3>, double,
double>(
point(Value), normal, Min_Curvature, Max_Curvature);
643 BRepAdaptor_Surface surf(face);
644 const double u0 = surf.FirstUParameter();
645 const double u1 = surf.LastUParameter();
646 const double v0 = surf.FirstVParameter();
647 const double v1 = surf.LastVParameter();
649 std::vector<CellData<2> > cells;
650 std::vector<Point<3> > vertices;
653 vertices.push_back(
point(surf.Value(u0,v0)));
654 vertices.push_back(
point(surf.Value(u1,v0)));
655 vertices.push_back(
point(surf.Value(u0,v1)));
656 vertices.push_back(
point(surf.Value(u1,v1)));
659 for (
unsigned int i=0; i<4; ++i)
662 cells.push_back(cell);
668 DEAL_II_NAMESPACE_CLOSE
TopoDS_Shape read_IGES(const std::string &filename, const double scale_factor=1e-3)
std_cxx11::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)
Point< 3 > line_intersection(const TopoDS_Shape &in_shape, const Point< 3 > &origin, const Tensor< 1, 3 > &direction, const double tolerance=1e-7)
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)
TopoDS_Shape read_STEP(const std::string &filename, const double scale_factor=1e-3)
#define AssertThrow(cond, exc)
numbers::NumberTraits< Number >::real_type norm() const
Point< 3 > push_forward(const TopoDS_Shape &in_shape, const double u, const double v)
void write_IGES(const TopoDS_Shape &shape, const std::string &filename)
TopoDS_Edge interpolation_curve(std::vector< Point< 3 > > &curve_points, const Tensor< 1, 3 > &direction=Tensor< 1, 3 >(), const bool closed=false, const double tolerance=1e-7)
std_cxx11::tuple< unsigned int, unsigned int, unsigned int > count_elements(const TopoDS_Shape &shape)
static ::ExceptionBase & ExcMessage(std::string arg1)
void write_STEP(const TopoDS_Shape &shape, const std::string &filename)
#define Assert(cond, exc)
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
std_cxx11::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)
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
std_cxx11::tuple< Point< 3 >, TopoDS_Shape, double, double > project_point_and_pull_back(const TopoDS_Shape &in_shape, const Point< 3 > &origin, const double tolerance=1e-7)
Point< 3 > closest_point(const TopoDS_Shape &in_shape, const Point< 3 > &origin, const double tolerance=1e-7)
double get_shape_tolerance(const TopoDS_Shape &shape)
TopoDS_Edge join_edges(const TopoDS_Shape &in_shape, 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)
static ::ExceptionBase & ExcUnsupportedShape()
Point< 3 > point(const gp_Pnt &p)
bool point_compare(const Point< 3 > &p1, const Point< 3 > &p2, const Tensor< 1, 3 > &direction=Tensor< 1, 3 >(), const double tolerance=1e-10)
void create_triangulation(const TopoDS_Face &face, Triangulation< 2, 3 > &tria)
unsigned int vertices[GeometryInfo< structdim >::vertices_per_cell]