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