WaveBEM: Unsteady Nonlinear Potential Flow Solver for Ship-Wave Interaction.
occ_utilities.h
Go to the documentation of this file.
1 #ifndef occ_utilities_h
2 #define occ_utilities_h
3 
4 #include <string>
5 #include <TopoDS_Shape.hxx>
6 
7 #include <Geom_Plane.hxx>
8 #include <Geom_Curve.hxx>
9 #include <gp_Pnt.hxx>
10 
11 #include <deal.II/base/point.h>
12 
13 
18 namespace OpenCascade
19 {
21  {
22  bool operator() (std::pair<unsigned int, double> a,std::pair<unsigned int, double> b) const
23  {
24  return a.second<b.second;
25  }
26  };
27 
29  // their content into openCascade
30  // topological entities. The option
31  // scale_factor is used to
32  // compensate for different units
33  // being used in the IGES files and
34  // in the target application.
35  TopoDS_Shape read_IGES(std::string filename, double scale_factor=1e-3);
36 
38  // given topological shape with
39  // the given plane. The returned
40  // topological shape will contain
41  // as few bsplines as possible.
42  void intersect_plane(const TopoDS_Shape &in_shape,
43  TopoDS_Shape &out_shape,
44  const double c_x,
45  const double c_y,
46  const double c_z,
47  const double c,
48  const double tolerance=1e-7);
49 
50  TopoDS_Shape extract_xz_edges(const TopoDS_Shape &in_shape,
51  const double tolerance=1e-7,
52  const unsigned int max_num_edges=20);
53 
54  // identifies transom edge, which is the
55  // edge with the highest mean x coordinate
56  // fails when transom is composed by more than
57  // one face/edge, so this method needs a trim
58  TopoDS_Shape extract_transom_edges(const TopoDS_Shape &in_shape,
59  const unsigned int num_transom_edges,
60  const double tolerance=1e-7);
61 
62  // creates a 3D smooth curve (bspline) passing
63  // through the points in the
64  // assigned vector. The points are
65  // reordered internally according
66  // to their scalar product with the
67  // direction, if direction is
68  // different from zero, otherwise
69  // they are used as passed.
70  Handle_Geom_Curve interpolation_curve_points_sort(std::vector<dealii::Point<3> > &curve_points,
71  dealii::Point<3> direction=dealii::Point<3>());
72 
73  Handle_Geom_Curve interpolation_curve(std::vector<dealii::Point<3> > &curve_points);
75  // Deal.II point.
76  inline dealii::Point<3> Pnt(const gp_Pnt &p)
77  {
78  dealii::Point<3> P(p.X(), p.Y(), p.Z());
79  return P;
80  }
81 
83  // OpenCascade point.
84  inline gp_Pnt Pnt(const dealii::Point<3> &p)
85  {
86  gp_Pnt P(p(0), p(1), p(2));
87  return P;
88  }
89 
90 
91  inline bool point_compare(const dealii::Point<3> &p1, const dealii::Point<3> &p2,
92  const dealii::Point<3> &direction)
93  {
94  return (p1*direction < p2*direction);
95  }
96 
98 }
99 
100 #endif
Handle_Geom_Curve interpolation_curve(std::vector< Point< 3 > > &curve_points)
dealii::Point< 3 > Pnt(const gp_Pnt &p)
Convert OpenCascade point into.
Definition: occ_utilities.h:76
void feature_edges_detection()
TopoDS_Shape extract_transom_edges(const TopoDS_Shape &in_shape, const unsigned int num_transom_edges, const double tolerance)
We collect in this namespace all utilities which operate on OpenCascade entities which don't need cla...
TopoDS_Shape read_IGES(string filename, double scale_factor)
Read IGES files and translate.
void intersect_plane(const TopoDS_Shape &in_shape, TopoDS_Shape &out_shape, const double c_x, const double c_y, const double c_z, const double c, const double tolerance)
Perform the intersection of the.
TopoDS_Shape extract_xz_edges(const TopoDS_Shape &in_shape, const double tolerance, const unsigned int max_num_edges)
bool point_compare(const dealii::Point< 3 > &p1, const dealii::Point< 3 > &p2, const dealii::Point< 3 > &direction)
Definition: occ_utilities.h:91
bool operator()(std::pair< unsigned int, double > a, std::pair< unsigned int, double > b) const
Definition: occ_utilities.h:22
Handle_Geom_Curve interpolation_curve_points_sort(std::vector< Point< 3 > > &curve_points, Point< 3 > direction)