|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
17 #ifndef dealii_occ_utilities_h
18 # define dealii_occ_utilities_h
22 # ifdef DEAL_II_WITH_OPENCASCADE
33 # define HAVE_CONFIG_H
34 # include <IFSelect_ReturnStatus.hxx>
35 # include <TopoDS_CompSolid.hxx>
36 # include <TopoDS_Compound.hxx>
37 # include <TopoDS_Edge.hxx>
38 # include <TopoDS_Face.hxx>
39 # include <TopoDS_Shape.hxx>
40 # include <TopoDS_Shell.hxx>
41 # include <TopoDS_Solid.hxx>
42 # include <TopoDS_Vertex.hxx>
43 # include <TopoDS_Wire.hxx>
44 # include <gp_Pnt.hxx>
112 std::tuple<unsigned int, unsigned int, unsigned int>
123 read_IGES(
const std::string &filename,
const double scale_factor = 1
e-3);
129 write_IGES(
const TopoDS_Shape &shape,
const std::string &filename);
138 read_STL(
const std::string &filename);
157 const std::string & filename,
158 const double deflection,
159 const bool sew_different_faces =
false,
160 const double sewer_tolerance = 1
e-6,
161 const bool is_relative =
false,
162 const double angular_deflection = 0.5,
163 const bool in_parallel =
false);
174 read_STEP(
const std::string &filename,
const double scale_factor = 1
e-3);
181 write_STEP(
const TopoDS_Shape &shape,
const std::string &filename);
210 const double tolerance = 1
e-7);
220 join_edges(
const TopoDS_Shape &in_shape,
const double tolerance = 1
e-7);
243 const bool closed =
false,
244 const double tolerance = 1
e-7);
253 std::vector<TopoDS_Face> & faces,
254 std::vector<TopoDS_Edge> & edges,
255 std::vector<TopoDS_Vertex> &
vertices);
269 template <
int spacedim>
295 template <
int spacedim>
296 std::vector<TopoDS_Edge>
309 std::vector<TopoDS_Compound> & compounds,
310 std::vector<TopoDS_CompSolid> &compsolids,
311 std::vector<TopoDS_Solid> & solids,
312 std::vector<TopoDS_Shell> & shells,
313 std::vector<TopoDS_Wire> & wires);
330 std::tuple<Point<dim>, TopoDS_Shape,
double,
double>
333 const double tolerance = 1
e-7);
345 const double tolerance = 1
e-7);
357 push_forward(
const TopoDS_Shape &in_shape,
const double u,
const double v);
369 const double tolerance = 1
e-7);
382 const double tolerance = 1
e-7);
397 const double tolerance = 1
e-7);
407 template <
int spacedim>
409 point(
const gp_Pnt &p,
const double tolerance = 1
e-10);
415 template <
int spacedim>
431 const double tolerance = 1
e-10);
441 <<
"The point [ " << arg1 <<
" ] is not on the manifold.");
450 <<
"Projection of point [ " << arg1 <<
" ] failed.");
456 IFSelect_ReturnStatus,
457 <<
"An OpenCASCADE routine failed with return status "
474 # endif // DEAL_II_WITH_OPENCASCADE
476 #endif // dealii_occ_utilities_h
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)
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)
static ::ExceptionBase & ExcProjectionFailed(Point< dim > arg1)
std::vector< TopoDS_Edge > create_curves_from_triangulation_boundary(const Triangulation< 2, spacedim > &triangulation, const Mapping< 2, spacedim > &mapping=StaticMappingQ1< 2, spacedim >::mapping)
static ::ExceptionBase & ExcOCCError(IFSelect_ReturnStatus arg1)
TopoDS_Edge join_edges(const TopoDS_Shape &in_shape, const double tolerance=1e-7)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
TopoDS_Shape read_STEP(const std::string &filename, const double scale_factor=1e-3)
void write_STEP(const TopoDS_Shape &shape, const std::string &filename)
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< dim > closest_point(const TopoDS_Shape &in_shape, const Point< dim > &origin, 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)
double get_shape_tolerance(const TopoDS_Shape &shape)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Point< dim > line_intersection(const TopoDS_Shape &in_shape, const Point< dim > &origin, const Tensor< 1, dim > &direction, const double tolerance=1e-7)
TopoDS_Shape read_STL(const std::string &filename)
TopoDS_Shape read_IGES(const std::string &filename, const double scale_factor=1e-3)
Abstract base class for mapping classes.
Point< dim > push_forward(const TopoDS_Shape &in_shape, const double u, const double v)
static ::ExceptionBase & ExcPointNotOnManifold(Point< dim > arg1)
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)
#define DEAL_II_NAMESPACE_OPEN
void create_triangulation(const TopoDS_Face &face, Triangulation< 2, spacedim > &tria)
#define DeclException1(Exception1, type1, outsequence)
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)
static ::ExceptionBase & ExcEdgeIsDegenerate()
static ::ExceptionBase & ExcUnsupportedShape()
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 extract_geometrical_shapes(const TopoDS_Shape &shape, std::vector< TopoDS_Face > &faces, std::vector< TopoDS_Edge > &edges, std::vector< TopoDS_Vertex > &vertices)
#define DeclException0(Exception0)
std::tuple< unsigned int, unsigned int, unsigned int > count_elements(const TopoDS_Shape &shape)
void write_IGES(const TopoDS_Shape &shape, const std::string &filename)
const typename ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
#define DEAL_II_NAMESPACE_CLOSE
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)