Reference documentation for deal.II version 8.4.2
utilities.cc
1 #include <deal.II/opencascade/utilities.h>
2 
3 #ifdef DEAL_II_WITH_OPENCASCADE
4 
5 #include <deal.II/base/point.h>
6 #include <deal.II/base/utilities.h>
7 #include <deal.II/base/exceptions.h>
8 
9 #include <boost/bind.hpp>
10 
11 #include <cstdio>
12 #include <iostream>
13 #include <set>
14 
15 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
16 
17 #include <IGESControl_Controller.hxx>
18 #include <IGESControl_Reader.hxx>
19 #include <IGESControl_Writer.hxx>
20 
21 #include <STEPControl_Controller.hxx>
22 #include <STEPControl_Reader.hxx>
23 #include <STEPControl_Writer.hxx>
24 
25 #include <TopoDS.hxx>
26 #include <TopoDS_Shape.hxx>
27 #include <TopoDS_Face.hxx>
28 #include <TopoDS_Edge.hxx>
29 #include <TopExp_Explorer.hxx>
30 
31 #include <Handle_Standard_Transient.hxx>
32 
33 #include <TColStd_SequenceOfTransient.hxx>
34 #include <TColStd_HSequenceOfTransient.hxx>
35 #include <TColgp_HArray1OfPnt.hxx>
36 
37 #include <gp_Pnt.hxx>
38 #include <gp_Lin.hxx>
39 #include <gp_Vec.hxx>
40 #include <GeomAPI_ProjectPointOnSurf.hxx>
41 #include <GeomAPI_ProjectPointOnCurve.hxx>
42 #include <IntCurvesFace_ShapeIntersector.hxx>
43 
44 #include <BRepTools.hxx>
45 #include <BRepAdaptor_Curve.hxx>
46 #include <BRepAdaptor_HCurve.hxx>
47 #include <BRepAdaptor_HCompCurve.hxx>
48 #include <BRepAdaptor_Surface.hxx>
49 #include <BRep_Builder.hxx>
50 #include <BRepBuilderAPI_Transform.hxx>
51 #include <BRepBuilderAPI_MakeEdge.hxx>
52 #include <BRepAlgo_Section.hxx>
53 
54 #include <Geom_Plane.hxx>
55 #include <Geom_BoundedCurve.hxx>
56 #include <GeomAPI_Interpolate.hxx>
57 #include <GeomConvert_CompCurveToBSplineCurve.hxx>
58 #include <GeomLProp_SLProps.hxx>
59 
60 #include <GCPnts_AbscissaPoint.hxx>
61 #include <ShapeAnalysis_Surface.hxx>
62 
64 
65 #include <vector>
66 #include <algorithm>
67 
68 DEAL_II_NAMESPACE_OPEN
69 
70 namespace OpenCASCADE
71 {
72  std_cxx11::tuple<unsigned int, unsigned int, unsigned int>
73  count_elements(const TopoDS_Shape &shape)
74  {
75  TopExp_Explorer exp;
76  unsigned int n_faces=0, n_edges=0, n_vertices=0;
77  for (exp.Init(shape, TopAbs_FACE);
78  exp.More(); exp.Next(), ++n_faces)
79  {}
80  for (exp.Init(shape, TopAbs_EDGE);
81  exp.More(); exp.Next(), ++n_edges)
82  {}
83  for (exp.Init(shape, TopAbs_VERTEX);
84  exp.More(); exp.Next(), ++n_vertices)
85  {}
86  return std_cxx11::tuple<unsigned int, unsigned int, unsigned int>(n_faces, n_edges, n_vertices);
87  }
88 
89  void extract_geometrical_shapes(const TopoDS_Shape &shape,
90  std::vector<TopoDS_Face> &faces,
91  std::vector<TopoDS_Edge> &edges,
92  std::vector<TopoDS_Vertex> &vertices)
93  {
94  faces.resize(0);
95  edges.resize(0);
96  vertices.resize(0);
97 
98  TopExp_Explorer exp;
99  for (exp.Init(shape, TopAbs_FACE); exp.More(); exp.Next())
100  {
101  faces.push_back(TopoDS::Face(exp.Current()));
102  }
103  for (exp.Init(shape, TopAbs_EDGE); exp.More(); exp.Next())
104  {
105  edges.push_back(TopoDS::Edge(exp.Current()));
106  }
107  for (exp.Init(shape, TopAbs_VERTEX); exp.More(); exp.Next())
108  {
109  vertices.push_back(TopoDS::Vertex(exp.Current()));
110  }
111  }
112 
113 
114  void extract_compound_shapes(const TopoDS_Shape &shape,
115  std::vector<TopoDS_Compound> &compounds,
116  std::vector<TopoDS_CompSolid> &compsolids,
117  std::vector<TopoDS_Solid> &solids,
118  std::vector<TopoDS_Shell> &shells,
119  std::vector<TopoDS_Wire> &wires)
120  {
121  compounds.resize(0);
122  compsolids.resize(0);
123  solids.resize(0);
124  shells.resize(0);
125  wires.resize(0);
126 
127  TopExp_Explorer exp;
128  for (exp.Init(shape, TopAbs_COMPOUND); exp.More(); exp.Next())
129  {
130  compounds.push_back(TopoDS::Compound(exp.Current()));
131  }
132  for (exp.Init(shape, TopAbs_COMPSOLID); exp.More(); exp.Next())
133  {
134  compsolids.push_back(TopoDS::CompSolid(exp.Current()));
135  }
136  for (exp.Init(shape, TopAbs_SOLID); exp.More(); exp.Next())
137  {
138  solids.push_back(TopoDS::Solid(exp.Current()));
139  }
140  for (exp.Init(shape, TopAbs_SHELL); exp.More(); exp.Next())
141  {
142  shells.push_back(TopoDS::Shell(exp.Current()));
143  }
144  for (exp.Init(shape, TopAbs_WIRE); exp.More(); exp.Next())
145  {
146  wires.push_back(TopoDS::Wire(exp.Current()));
147  }
148  }
149 
150  gp_Pnt point(const Point<3> &p)
151  {
152  return gp_Pnt(p(0), p(1), p(2));
153  }
154 
155 
156  Point<3> point(const gp_Pnt &p)
157  {
158  return Point<3>(p.X(), p.Y(), p.Z());
159  }
160 
161  bool point_compare(const Point<3> &p1, const Point<3> &p2,
162  const Tensor<1,3> &direction,
163  const double tolerance)
164  {
165  const double rel_tol=std::max(tolerance, std::max(p1.norm(), p2.norm())*tolerance);
166  if (direction.norm() > 0.0)
167  return (p1*direction < p2*direction-rel_tol);
168  else
169  for (int d=2; d>=0; --d)
170  if (p1[d] < p2[d]-rel_tol)
171  return true;
172  else if (p2[d] < p1[d]-rel_tol)
173  return false;
174 
175  // If we got here, for all d, none of the conditions above was
176  // satisfied. The two points are equal up to tolerance
177  return false;
178  }
179 
180 
181  TopoDS_Shape read_IGES(const std::string &filename,
182  const double scale_factor)
183  {
184  IGESControl_Reader reader;
185  IFSelect_ReturnStatus stat;
186  stat = reader.ReadFile(filename.c_str());
187  AssertThrow(stat == IFSelect_RetDone,
188  ExcMessage("Error in reading file!"));
189 
190  Standard_Boolean failsonly = Standard_False;
191  IFSelect_PrintCount mode = IFSelect_ItemsByEntity;
192  reader.PrintCheckLoad (failsonly, mode);
193 
194  Standard_Integer nRoots = reader.TransferRoots();
195  //selects all IGES entities (including non visible ones) in the
196  //file and puts them into a list called MyList,
197 
198  AssertThrow(nRoots > 0,
199  ExcMessage("Read nothing from file."));
200 
201  // Handle IGES Scale here.
202  gp_Pnt Origin;
203  gp_Trsf scale;
204  scale.SetScale (Origin, scale_factor);
205 
206  TopoDS_Shape sh = reader.OneShape();
207  BRepBuilderAPI_Transform trans(sh, scale);
208 
209  return trans.Shape(); // this is the actual translation
210  }
211 
212  void write_IGES(const TopoDS_Shape &shape,
213  const std::string &filename)
214  {
215  IGESControl_Controller::Init();
216  IGESControl_Writer ICW ("MM", 0);
217  Standard_Boolean ok = ICW.AddShape (shape);
218  AssertThrow(ok, ExcMessage("Failed to add shape to IGES controller."));
219  ICW.ComputeModel();
220  Standard_Boolean OK = ICW.Write (filename.c_str());
221  AssertThrow(OK, ExcMessage("Failed to write IGES file."));
222  }
223 
224  TopoDS_Shape read_STEP(const std::string &filename,
225  const double scale_factor)
226  {
227  STEPControl_Reader reader;
228  IFSelect_ReturnStatus stat;
229  stat = reader.ReadFile(filename.c_str());
230  AssertThrow(stat == IFSelect_RetDone,
231  ExcMessage("Error in reading file!"));
232 
233  Standard_Boolean failsonly = Standard_False;
234  IFSelect_PrintCount mode = IFSelect_ItemsByEntity;
235  reader.PrintCheckLoad (failsonly, mode);
236 
237  Standard_Integer nRoots = reader.TransferRoots();
238  //selects all IGES entities (including non visible ones) in the
239  //file and puts them into a list called MyList,
240 
241  AssertThrow(nRoots > 0,
242  ExcMessage("Read nothing from file."));
243 
244  // Handle STEP Scale here.
245  gp_Pnt Origin;
246  gp_Trsf scale;
247  scale.SetScale (Origin, scale_factor);
248 
249  TopoDS_Shape sh = reader.OneShape();
250  BRepBuilderAPI_Transform trans(sh, scale);
251 
252  return trans.Shape(); // this is the actual translation
253  }
254 
255  void write_STEP(const TopoDS_Shape &shape,
256  const std::string &filename)
257  {
258  STEPControl_Controller::Init();
259  STEPControl_Writer SCW;
260  IFSelect_ReturnStatus status;
261  status = SCW.Transfer(shape, STEPControl_AsIs);
262  AssertThrow(status == IFSelect_RetDone, ExcMessage("Failed to add shape to STEP controller."));
263 
264  status = SCW.Write(filename.c_str());
265 
266  AssertThrow(status == IFSelect_RetDone, ExcMessage("Failed to write translated shape to STEP file."));
267  }
268 
269  double get_shape_tolerance(const TopoDS_Shape &shape)
270  {
271  double tolerance = 0.0;
272 
273  std::vector<TopoDS_Face> faces;
274  std::vector<TopoDS_Edge> edges;
275  std::vector<TopoDS_Vertex> vertices;
276 
278  faces,
279  edges,
280  vertices);
281 
282  for (unsigned int i=0; i<vertices.size(); ++i)
283  tolerance = fmax(tolerance,BRep_Tool::Tolerance(vertices[i]));
284 
285  for (unsigned int i=0; i<edges.size(); ++i)
286  tolerance = fmax(tolerance,BRep_Tool::Tolerance(edges[i]));
287 
288  for (unsigned int i=0; i<faces.size(); ++i)
289  tolerance = fmax(tolerance,BRep_Tool::Tolerance(faces[i]));
290 
291 
292  return tolerance;
293  }
294 
295  TopoDS_Shape intersect_plane(const TopoDS_Shape &in_shape,
296  const double c_x,
297  const double c_y,
298  const double c_z,
299  const double c,
300  const double /*tolerance*/)
301  {
302  Handle(Geom_Plane) plane = new Geom_Plane(c_x,c_y,c_z,c);
303  BRepAlgo_Section section(in_shape, plane);
304  TopoDS_Shape edges = section.Shape();
305  return edges;
306  }
307 
308  TopoDS_Edge join_edges(const TopoDS_Shape &in_shape,
309  const double tolerance)
310  {
311  TopoDS_Edge out_shape;
312  TopoDS_Shape edges = in_shape;
313  std::vector<Handle_Geom_BoundedCurve> intersections;
314  TopLoc_Location L;
315  Standard_Real First;
316  Standard_Real Last;
317  gp_Pnt PIn(0.0,0.0,0.0);
318  gp_Pnt PFin(0.0,0.0,0.0);
319  gp_Pnt PMid(0.0,0.0,0.0);
320  TopExp_Explorer edgeExplorer(edges , TopAbs_EDGE);
321  TopoDS_Edge edge;
322  while (edgeExplorer.More())
323  {
324  edge = TopoDS::Edge(edgeExplorer.Current());
325  Handle(Geom_Curve) curve = BRep_Tool::Curve(edge,L,First,Last);
326  intersections.push_back(Handle(Geom_BoundedCurve)::DownCast(curve));
327  edgeExplorer.Next();
328  }
329 
330  // Now we build a single bspline out of all the geometrical
331  // curves, in Lexycographical order
332  unsigned int numIntersEdges = intersections.size();
333  Assert(numIntersEdges>0, ExcMessage("No curves to process!"));
334 
335  GeomConvert_CompCurveToBSplineCurve convert_bspline(intersections[0]);
336 
337  bool check = false, one_added = true, one_failed=true;
338  std::vector<bool> added(numIntersEdges, false);
339  added[0] = true;
340  while (one_added == true)
341  {
342  one_added = false;
343  one_failed = false;
344  for (unsigned int i=1; i<numIntersEdges; ++i)
345  if (added[i] == false)
346  {
347  Handle(Geom_Curve) curve = intersections[i];
348  Handle(Geom_BoundedCurve) bcurve = Handle(Geom_BoundedCurve)::DownCast(curve);
349  check = convert_bspline.Add(bcurve,tolerance,0,1,0);
350  if (check == false) // If we failed, try again with the reversed curve
351  {
352  curve->Reverse();
353  Handle(Geom_BoundedCurve) bcurve = Handle(Geom_BoundedCurve)::DownCast(curve);
354  check = convert_bspline.Add(bcurve,tolerance,0,1,0);
355  }
356  one_failed = one_failed || (check == false);
357  one_added = one_added || (check == true);
358  added[i] = check;
359  }
360  }
361 
362  Assert(one_failed == false,
363  ExcMessage("Joining some of the Edges failed."));
364 
365  Handle(Geom_Curve) bspline = convert_bspline.BSplineCurve();
366 
367  out_shape = BRepBuilderAPI_MakeEdge(bspline);
368  return out_shape;
369  }
370 
371 
372  Point<3> line_intersection(const TopoDS_Shape &in_shape,
373  const Point<3> &origin,
374  const Tensor<1,3> &direction,
375  const double tolerance)
376  {
377  // translating original Point<dim> to gp point
378 
379  gp_Pnt P0 = point(origin);
380  gp_Ax1 gpaxis(P0, gp_Dir(direction[0], direction[1], direction[2]));
381  gp_Lin line(gpaxis);
382 
383  // destination point
384  gp_Pnt Pproj(0.0,0.0,0.0);
385 
386  // we prepare now the surface for the projection we get the whole
387  // shape from the iges model
388  IntCurvesFace_ShapeIntersector Inters;
389  Inters.Load(in_shape,tolerance);
390 
391  // Keep in mind: PerformNearest sounds pretty but DOESN'T WORK!!!
392  // The closest point must be found by hand
393  Inters.Perform(line,-RealLast(),+RealLast());
394  Assert(Inters.IsDone(), ExcMessage("Could not project point."));
395 
396  double minDistance = 1e7;
397  double distance;
398  Point<3> result;
399  for (int i=0; i<Inters.NbPnt(); ++i)
400  {
401  distance = point(origin).Distance(Inters.Pnt(i+1));
402  //cout<<"Point "<<i<<": "<<point(Inters.Pnt(i+1))<<" distance: "<<distance<<endl;
403  if (distance < minDistance)
404  {
405  minDistance = distance;
406  result = point(Inters.Pnt(i+1));
407  }
408  }
409 
410  return result;
411  }
412 
413  TopoDS_Edge interpolation_curve(std::vector<Point<3> > &curve_points,
414  const Tensor<1,3> &direction,
415  const bool closed,
416  const double tolerance)
417  {
418 
419  unsigned int n_vertices = curve_points.size();
420 
421  if (direction*direction > 0)
422  {
423  std::sort(curve_points.begin(), curve_points.end(),
424  boost::bind(&OpenCASCADE::point_compare, _1, _2, direction, tolerance));
425  }
426 
427  // set up array of vertices
428  Handle(TColgp_HArray1OfPnt) vertices = new TColgp_HArray1OfPnt(1,n_vertices);
429  for (unsigned int vertex=0; vertex<n_vertices; ++vertex)
430  {
431  vertices->SetValue(vertex+1,point(curve_points[vertex]));
432  }
433 
434 
435  GeomAPI_Interpolate bspline_generator(vertices, closed, tolerance);
436  bspline_generator.Perform();
437  Assert( (bspline_generator.IsDone()), ExcMessage("Interpolated bspline generation failed"));
438 
439  Handle(Geom_BSplineCurve) bspline = bspline_generator.Curve();
440  TopoDS_Edge out_shape = BRepBuilderAPI_MakeEdge(bspline);
441  return out_shape;
442  }
443 
444  std_cxx11::tuple<Point<3>, TopoDS_Shape, double, double>
445  project_point_and_pull_back(const TopoDS_Shape &in_shape,
446  const Point<3> &origin,
447  const double tolerance)
448  {
449  TopExp_Explorer exp;
450  gp_Pnt Pproj = point(origin);
451 
452  double minDistance = 1e7;
453  gp_Pnt tmp_proj(0.0,0.0,0.0);
454 
455  unsigned int counter = 0;
456  unsigned int face_counter = 0;
457 
458  TopoDS_Shape out_shape;
459  double u=0;
460  double v=0;
461 
462  for (exp.Init(in_shape, TopAbs_FACE); exp.More(); exp.Next())
463  {
464  TopoDS_Face face = TopoDS::Face(exp.Current());
465 
466  // the projection function needs a surface, so we obtain the
467  // surface upon which the face is defined
468  Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
469 
470  ShapeAnalysis_Surface projector(SurfToProj);
471  gp_Pnt2d proj_params = projector.ValueOfUV(point(origin), tolerance);
472 
473  SurfToProj->D0(proj_params.X(),proj_params.Y(),tmp_proj);
474 
475  double distance = point(tmp_proj).distance(origin);
476  if (distance < minDistance)
477  {
478  minDistance = distance;
479  Pproj = tmp_proj;
480  out_shape = face;
481  u=proj_params.X();
482  v=proj_params.Y();
483  ++counter;
484  }
485  ++face_counter;
486  }
487 
488  // face counter tells us if the shape contained faces: if it does, there is no need
489  // to loop on edges. Even if the closest point lies on the boundary of a parametric surface,
490  // we need in fact to retain the face and both u and v, if we want to use this method to
491  // retrieve the surface normal
492  if (face_counter==0)
493  for (exp.Init(in_shape, TopAbs_EDGE); exp.More(); exp.Next())
494  {
495  TopoDS_Edge edge = TopoDS::Edge(exp.Current());
496  if (!BRep_Tool::Degenerated(edge))
497  {
498  TopLoc_Location L;
499  Standard_Real First;
500  Standard_Real Last;
501 
502  // the projection function needs a Curve, so we obtain the
503  // curve upon which the edge is defined
504  Handle(Geom_Curve) CurveToProj = BRep_Tool::Curve(edge,L,First,Last);
505 
506  GeomAPI_ProjectPointOnCurve Proj(point(origin),CurveToProj);
507  unsigned int num_proj_points = Proj.NbPoints();
508  if ((num_proj_points > 0) && (Proj.LowerDistance() < minDistance))
509  {
510  minDistance = Proj.LowerDistance();
511  Pproj = Proj.NearestPoint();
512  out_shape = edge;
513  u=Proj.LowerDistanceParameter();
514  ++counter;
515  }
516  }
517  }
518 
519  Assert(counter > 0, ExcMessage("Could not find projection points."));
520  return std_cxx11::tuple<Point<3>, TopoDS_Shape, double, double>
521  (point(Pproj),out_shape, u, v);
522  }
523 
524 
525  Point<3> closest_point(const TopoDS_Shape &in_shape,
526  const Point<3> &origin,
527  const double tolerance)
528  {
529  std_cxx11::tuple<Point<3>, TopoDS_Shape, double, double>
530  ref = project_point_and_pull_back(in_shape, origin, tolerance);
531  return std_cxx11::get<0>(ref);
532  }
533 
534  std_cxx11::tuple<Point<3>, Point<3>, double>
535  closest_point_and_differential_forms(const TopoDS_Shape &in_shape,
536  const Point<3> &origin,
537  const double tolerance)
538 
539  {
540  std_cxx11::tuple<Point<3>, TopoDS_Shape, double, double>
541  shape_and_params = project_point_and_pull_back(in_shape,
542  origin,
543  tolerance);
544 
545  TopoDS_Shape &out_shape = std_cxx11::get<1>(shape_and_params);
546  double &u = std_cxx11::get<2>(shape_and_params);
547  double &v = std_cxx11::get<3>(shape_and_params);
548 
549  // just a check here: the number of faces in out_shape must be 1, otherwise
550  // something is wrong
551  std_cxx11::tuple<unsigned int, unsigned int, unsigned int> numbers =
552  count_elements(out_shape);
553  (void)numbers;
554 
555  Assert(std_cxx11::get<0>(numbers) > 0,
556  ExcMessage("Could not find normal: the shape containing the closest point has 0 faces."));
557  Assert(std_cxx11::get<0>(numbers) < 2,
558  ExcMessage("Could not find normal: the shape containing the closest point has more than 1 face."));
559 
560 
561  TopExp_Explorer exp;
562  exp.Init(out_shape, TopAbs_FACE);
563  TopoDS_Face face = TopoDS::Face(exp.Current());
564  return push_forward_and_differential_forms(face, u, v, tolerance);
565  }
566 
567  Point<3> push_forward(const TopoDS_Shape &in_shape,
568  const double u,
569  const double v)
570  {
571  switch (in_shape.ShapeType())
572  {
573  case TopAbs_FACE:
574  {
575  BRepAdaptor_Surface surf(TopoDS::Face(in_shape));
576  return point(surf.Value(u,v));
577  }
578  case TopAbs_EDGE:
579  {
580  BRepAdaptor_Curve curve(TopoDS::Edge(in_shape));
581  return point(curve.Value(u));
582  }
583  default:
584  Assert(false, ExcUnsupportedShape());
585  }
586  return Point<3>();
587  }
588 
589  std_cxx11::tuple<Point<3>, Point<3>, double >
590  push_forward_and_differential_forms(const TopoDS_Face &face,
591  const double u,
592  const double v,
593  const double /*tolerance*/)
594  {
595  Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
596  GeomLProp_SLProps props(SurfToProj, u, v, 1, 1e-7);
597  gp_Pnt Value = props.Value();
598  Assert(props.IsNormalDefined(), ExcMessage("Normal is not well defined!"));
599  gp_Dir Normal = props.Normal();
600  Assert(props.IsCurvatureDefined(), ExcMessage("Curvature is not well defined!"));
601  Standard_Real Mean_Curvature = props.MeanCurvature();
602  Point<3> normal = Point<3>(Normal.X(),Normal.Y(),Normal.Z());
603  return std_cxx11::tuple<Point<3>, Point<3>, double>(point(Value), normal, Mean_Curvature);
604  }
605 
606 
607 
608  void create_triangulation(const TopoDS_Face &face,
609  Triangulation<2,3> &tria)
610  {
611  BRepAdaptor_Surface surf(face);
612  const double u0 = surf.FirstUParameter();
613  const double u1 = surf.LastUParameter();
614  const double v0 = surf.FirstVParameter();
615  const double v1 = surf.LastVParameter();
616 
617  std::vector<CellData<2> > cells;
618  std::vector<Point<3> > vertices;
619  SubCellData t;
620 
621  vertices.push_back(point(surf.Value(u0,v0)));
622  vertices.push_back(point(surf.Value(u1,v0)));
623  vertices.push_back(point(surf.Value(u0,v1)));
624  vertices.push_back(point(surf.Value(u1,v1)));
625 
626  CellData<2> cell;
627  for (unsigned int i=0; i<4; ++i)
628  cell.vertices[i] = i;
629 
630  cells.push_back(cell);
631  tria.create_triangulation(vertices, cells, t);
632  }
633 
634 } // end namespace
635 
636 DEAL_II_NAMESPACE_CLOSE
637 
638 #endif
std_cxx11::tuple< Point< 3 >, Point< 3 >, double > closest_point_and_differential_forms(const TopoDS_Shape &in_shape, const Point< 3 > &origin, const double tolerance=1e-7)
Definition: utilities.cc:535
TopoDS_Shape read_IGES(const std::string &filename, const double scale_factor=1e-3)
Definition: utilities.cc:181
Point< 3 > line_intersection(const TopoDS_Shape &in_shape, const Point< 3 > &origin, const Tensor< 1, 3 > &direction, const double tolerance=1e-7)
Definition: utilities.cc:372
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)
Definition: utilities.cc:295
TopoDS_Shape read_STEP(const std::string &filename, const double scale_factor=1e-3)
Definition: utilities.cc:224
::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
numbers::NumberTraits< Number >::real_type norm() const
Definition: tensor.h:1047
Point< 3 > push_forward(const TopoDS_Shape &in_shape, const double u, const double v)
Definition: utilities.cc:567
void write_IGES(const TopoDS_Shape &shape, const std::string &filename)
Definition: utilities.cc:212
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)
Definition: utilities.cc:413
std_cxx11::tuple< unsigned int, unsigned int, unsigned int > count_elements(const TopoDS_Shape &shape)
Definition: utilities.cc:73
std_cxx11::tuple< Point< 3 >, Point< 3 >, double > push_forward_and_differential_forms(const TopoDS_Face &face, const double u, const double v, const double tolerance=1e-7)
Definition: utilities.cc:590
void write_STEP(const TopoDS_Shape &shape, const std::string &filename)
Definition: utilities.cc:255
#define Assert(cond, exc)
Definition: exceptions.h:294
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
Definition: tria.cc:9130
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)
Definition: utilities.cc:114
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)
Definition: utilities.cc:445
Point< 3 > closest_point(const TopoDS_Shape &in_shape, const Point< 3 > &origin, const double tolerance=1e-7)
Definition: utilities.cc:525
double get_shape_tolerance(const TopoDS_Shape &shape)
Definition: utilities.cc:269
TopoDS_Edge join_edges(const TopoDS_Shape &in_shape, const double tolerance=1e-7)
Definition: utilities.cc:308
Definition: mpi.h:48
void extract_geometrical_shapes(const TopoDS_Shape &shape, std::vector< TopoDS_Face > &faces, std::vector< TopoDS_Edge > &edges, std::vector< TopoDS_Vertex > &vertices)
Definition: utilities.cc:89
Point< 3 > point(const gp_Pnt &p)
Definition: utilities.cc:156
bool point_compare(const Point< 3 > &p1, const Point< 3 > &p2, const Tensor< 1, 3 > &direction=Tensor< 1, 3 >(), const double tolerance=1e-10)
Definition: utilities.cc:161
void create_triangulation(const TopoDS_Face &face, Triangulation< 2, 3 > &tria)
Definition: utilities.cc:608
unsigned int vertices[GeometryInfo< structdim >::vertices_per_cell]
Definition: tria.h:110