Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
Classes | Functions
OpenCASCADE Namespace Reference

Classes

class  ArclengthProjectionLineManifold
 
class  DirectionalProjectionManifold
 
class  NormalProjectionManifold
 
class  NormalToMeshProjectionManifold
 
class  NURBSPatchManifold
 

Functions

std::tuple< unsigned int, unsigned int, unsigned intcount_elements (const TopoDS_Shape &shape)
 
TopoDS_Shape read_IGES (const std::string &filename, const double scale_factor=1e-3)
 
void write_IGES (const TopoDS_Shape &shape, const std::string &filename)
 
TopoDS_Shape read_STL (const std::string &filename)
 
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)
 
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)
 
double get_shape_tolerance (const TopoDS_Shape &shape)
 
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_Edge join_edges (const TopoDS_Shape &in_shape, const double tolerance=1e-7)
 
template<int dim>
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)
 
void extract_geometrical_shapes (const TopoDS_Shape &shape, std::vector< TopoDS_Face > &faces, std::vector< TopoDS_Edge > &edges, std::vector< TopoDS_Vertex > &vertices)
 
template<int spacedim>
void create_triangulation (const TopoDS_Face &face, Triangulation< 2, spacedim > &tria)
 
template<int spacedim>
std::vector< TopoDS_Edge > create_curves_from_triangulation_boundary (const Triangulation< 2, spacedim > &triangulation, const Mapping< 2, spacedim > &mapping=StaticMappingQ1< 2, spacedim >::mapping)
 
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)
 
template<int dim>
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)
 
template<int dim>
Point< dim > closest_point (const TopoDS_Shape &in_shape, const Point< dim > &origin, const double tolerance=1e-7)
 
template<int dim>
Point< dim > push_forward (const TopoDS_Shape &in_shape, const double u, const double v)
 
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)
 
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)
 
template<int dim>
Point< dim > line_intersection (const TopoDS_Shape &in_shape, const Point< dim > &origin, const Tensor< 1, dim > &direction, const double tolerance=1e-7)
 
template<int spacedim>
Point< spacedim > point (const gp_Pnt &p, const double tolerance=1e-10)
 
template<int spacedim>
gp_Pnt point (const Point< spacedim > &p)
 
template<int dim>
bool point_compare (const Point< dim > &p1, const Point< dim > &p2, const Tensor< 1, dim > &direction=Tensor< 1, dim >(), const double tolerance=1e-10)
 
template<int dim>
static ::ExceptionBaseExcPointNotOnManifold (Point< dim > arg1)
 
template<int dim>
static ::ExceptionBaseExcProjectionFailed (Point< dim > arg1)
 
static ::ExceptionBaseExcOCCError (IFSelect_ReturnStatus arg1)
 
static ::ExceptionBaseExcEdgeIsDegenerate ()
 
static ::ExceptionBaseExcUnsupportedShape ()
 

Detailed Description

We collect in this namespace all utilities which operate on OpenCASCADE entities. OpenCASCADE splits every object into a topological description and a geometrical entity. The basic topological description is a TopoDS_Shape. TopoDS_Shapes are light objects, and can be copied around. The closest deal.II analog is a TriaIterator.

The OpenCASCADE topology is designed with reference to the STEP standard ISO-10303-42. The structure is an oriented one-way graph, where parents refer to their children, and there are no back references. Abstract structure is implemented as C++ classes from the TopoDS package. A TopoDS_Shape is manipulated by value and contains 3 fields: location, orientation and a myTShape handle (of the TopoDS_TShape type). According to OpenCASCADE documentation, myTShape and Location are used to share data between various shapes to save memory. For example, an edge belonging to two faces has equal Locations and myTShape fields but different Orientations (Forward in context of one face and Reversed in one of the other).

Valid shapes include collection of other shapes, solids, faces, edges, vertices, etc.

Once a topological description is available, if a concrete geometrical object can be created, the BRep classes allow one to extract the actual geometrical information from a shape.

This is done by inheriting abstract topology classes from the TopoDS package by those implementing a boundary representation model (from the BRep package). Only 3 types of topological objects have geometric representations - vertex, edge, and face.

Every TopoDS_Shape can be queried to figure out what type of shape it is, and actual geometrical objects, like surfaces, curves or points, can be extracted using BRepTools.

In this namespace we provide readers and writers that read standard CAD files, and return a TopoDS_Shape, or that write a CAD file, given a TopoDS_Shape. Most of the functions in the OpenCASCADE namespace deal with TopoDS_Shapes of one type or another, and provide interfaces to common deal.II objects, like Triangulation, Manifold, and so on.

Notice that most of these tools are only useful when spacedim is equal to three, since OpenCASCADE only operates in three-dimensional mode. In some cases they can be used in two dimensions as well, and the third dimension will be set to zero.

If you wish to use these tools when the dimension of the space is two, then make sure your CAD files are actually flat and that all z coordinates are equal to zero, as otherwise you will get many exceptions.

Function Documentation

◆ count_elements()

std::tuple< unsigned int, unsigned int, unsigned int > OpenCASCADE::count_elements ( const TopoDS_Shape &  shape)

Count the subobjects of a shape. This function is useful to gather information about the TopoDS_Shape passed as argument. It returns the number of faces, edges and vertices (the only topological entities associated with actual geometries) which are contained in the given shape.

Definition at line 92 of file utilities.cc.

◆ read_IGES()

TopoDS_Shape OpenCASCADE::read_IGES ( const std::string &  filename,
const double  scale_factor = 1e-3 
)

Read IGES files and translate their content into openCascade topological entities. The option scale_factor is used to compensate for different units being used in the IGES files and in the target application. The standard unit for IGES files is millimeters. The return object is a TopoDS_Shape which contains all objects from the file.

Definition at line 239 of file utilities.cc.

◆ write_IGES()

void OpenCASCADE::write_IGES ( const TopoDS_Shape &  shape,
const std::string &  filename 
)

Write the given topological shape into an IGES file.

Definition at line 268 of file utilities.cc.

◆ read_STL()

TopoDS_Shape OpenCASCADE::read_STL ( const std::string &  filename)

Read STL files and translate their content into openCascade topological entities. The return object is a TopoDS_Shape which contains all objects from the file.

Definition at line 281 of file utilities.cc.

◆ write_STL()

void OpenCASCADE::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 
)

Write the given topological shape into an STL file. In order to do so the shape must contain a mesh structure, the function checks if all the faces of the shape have an attached mesh, if this is not the case it proceeds to mesh it automatically. We remark that the automatic mesh generation in OpenCASCADE takes care only of the geometrical resemblance between the shape and the mesh, to control the shape and regularity of the triangles you should use other meshing softwares. The two arguments deflection and angular_deflection select the accuracy of the created triangulation with respect to the original topological shape. The argument sew_different_faces gives the possibility to use a Sewer from OpenCASCADE to create a watertight closed STL using the argument sewer_tolerance. The argument is_relative specifies if distance are relative and in_parallel if the execution should be in parallel.

Definition at line 291 of file utilities.cc.

◆ read_STEP()

TopoDS_Shape OpenCASCADE::read_STEP ( const std::string &  filename,
const double  scale_factor = 1e-3 
)

Read STEP files and translate their content into openCascade topological entities. The option scale_factor is used to compensate for different units being used in the STEP files and in the target application. The standard unit for STEP files is millimeters. The return object is a TopoDS_Shape which contains all objects from the file.

Definition at line 356 of file utilities.cc.

◆ write_STEP()

void OpenCASCADE::write_STEP ( const TopoDS_Shape &  shape,
const std::string &  filename 
)

Write the given topological shape into an STEP file.

Definition at line 385 of file utilities.cc.

◆ get_shape_tolerance()

double OpenCASCADE::get_shape_tolerance ( const TopoDS_Shape &  shape)

This function returns the tolerance associated with the shape. Each CAD geometrical object is defined along with a tolerance, which indicates possible inaccuracy of its placement. For instance, the tolerance of a vertex indicates that it can be located in any point contained in a sphere centered in the nominal position and having radius tol. While carrying out an operation such as projecting a point onto a surface (which will in turn have its tolerance) we must keep in mind that the precision of the projection will be limited by the tolerance with which the surface is built. The tolerance is computed taking the maximum tolerance among the subshapes composing the shape.

Definition at line 401 of file utilities.cc.

◆ intersect_plane()

TopoDS_Shape OpenCASCADE::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 
)

Perform the intersection of the given topological shape with the plane \(c_x x + c_y y + c_z z +c = 0\). The returned topological shape will contain as few bsplines as possible. An exception is thrown if the intersection produces an empty shape.

Definition at line 425 of file utilities.cc.

◆ join_edges()

TopoDS_Edge OpenCASCADE::join_edges ( const TopoDS_Shape &  in_shape,
const double  tolerance = 1e-7 
)

Try to join all edges contained in the given TopoDS_Shape into a single TopoDS_Edge, containing as few BSPlines as possible. If the input shape contains faces, they will be ignored by this function. If the contained edges cannot be joined into a single one, i.e., they form disconnected curves, an exception will be thrown.

Definition at line 443 of file utilities.cc.

◆ interpolation_curve()

template<int dim>
TopoDS_Edge OpenCASCADE::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 
)

Creates a smooth BSpline curve passing through the points in the assigned vector, and store it in the returned TopoDS_Shape (which is of type TopoDS_Edge). The points are reordered internally according to their scalar product with the direction, if direction is different from zero, otherwise they are used as passed. Notice that this function changes the input points if required by the algorithm.

This class is used to interpolate a BsplineCurve passing through an array of points, with a C2 Continuity. If the optional parameter closed is set to true, then the curve will be C2 at all points except the first (where only C1 continuity will be given), and it will be a closed curve.

The curve is guaranteed to be at distance tolerance from the input points. If the algorithm fails in generating such a curve, an exception is thrown.

Definition at line 557 of file utilities.cc.

◆ extract_geometrical_shapes()

void OpenCASCADE::extract_geometrical_shapes ( const TopoDS_Shape &  shape,
std::vector< TopoDS_Face > &  faces,
std::vector< TopoDS_Edge > &  edges,
std::vector< TopoDS_Vertex > &  vertices 
)

Extract all subshapes from a TopoDS_Shape, and store the results into standard containers. If the shape does not contain a certain type of shape, the respective container will be empty.

Definition at line 108 of file utilities.cc.

◆ create_triangulation()

template<int spacedim>
void OpenCASCADE::create_triangulation ( const TopoDS_Face &  face,
Triangulation< 2, spacedim > &  tria 
)

Create a triangulation from a single face. This class extracts the first u and v parameter of the parametric surface making up this face, and creates a Triangulation<2,spacedim> containing a single coarse cell reflecting this face. If the surface is not a trimmed surface, the vertices of this cell will coincide with the TopoDS_Vertex vertices of the original TopoDS_Face. This, however, is often not the case, and the user should be careful on how this mesh is used.

If you call this function with a Triangulation<2,2>, make sure that the input face has all z coordinates set to zero, or you'll get an exception.

Definition at line 895 of file utilities.cc.

◆ create_curves_from_triangulation_boundary()

template<int spacedim>
std::vector< TopoDS_Edge > OpenCASCADE::create_curves_from_triangulation_boundary ( const Triangulation< 2, spacedim > &  triangulation,
const Mapping< 2, spacedim > &  mapping = StaticMappingQ1<2, spacedim>::mapping 
)

Given a Triangulation and an optional Mapping, create a vector of smooth curves that interpolate the connected parts of the boundary vertices of the Triangulation and return them as a vector of TopoDS_Edge objects.

This function constructs closed Bspline curve objects passing through all vertices of the boundary of the triangulation, with \(C^2\) Continuity on each vertex except the first, where only \(C^1\) continuity is guaranteed.

The returned curves are ordered with respect to the indices of the faces that make up the triangulation boundary, i.e., the first curve is the one extracted starting from the face with the lowest index, and so on.

Parameters
[in]triangulationInput triangulation
[in]mappingOptional input mapping
Returns
An std::vector of TopoDS_Edge objects representing the smooth interpolation of the boundary of the triangulation

Definition at line 600 of file utilities.cc.

◆ extract_compound_shapes()

void OpenCASCADE::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 
)

Extract all compound shapes from a TopoDS_Shape, and store the results into standard containers. If the shape does not contain a certain type of compound, the respective container will be empty.

Definition at line 134 of file utilities.cc.

◆ project_point_and_pull_back()

template<int dim>
std::tuple< Point< dim >, TopoDS_Shape, double, double > OpenCASCADE::project_point_and_pull_back ( const TopoDS_Shape &  in_shape,
const Point< dim > &  origin,
const double  tolerance = 1e-7 
)

Project the point origin on the topological shape given by in_shape, and returns the projected point, the subshape which contains the point and the parametric u and v coordinates of the point within the resulting shape. If the shape is not elementary, all its subshapes are iterated, faces first, then edges, and the returned shape is the closest one to the point origin. If the returned shape is an edge, then only the u coordinate is filled with sensible information, and the v coordinate is set to zero.

This function returns a tuple containing the projected point, the shape, the u coordinate and the v coordinate (which is different from zero only if the resulting shape is a face).

Definition at line 703 of file utilities.cc.

◆ closest_point()

template<int dim>
Point< dim > OpenCASCADE::closest_point ( const TopoDS_Shape &  in_shape,
const Point< dim > &  origin,
const double  tolerance = 1e-7 
)

Return the projection of the point origin on the topological shape given by in_shape. If the shape is not elementary, all its subshapes are iterated, faces first, then edges, and the returned point is the closest one to the in_shape, regardless of its type.

Definition at line 791 of file utilities.cc.

◆ push_forward()

template<int dim>
Point< dim > OpenCASCADE::push_forward ( const TopoDS_Shape &  in_shape,
const double  u,
const double  v 
)

Given an elementary shape in_shape and the reference coordinates within the shape, returns the corresponding point in real space. If the shape is a TopoDS_Edge, the v coordinate is ignored. Only edges or faces, as returned by the function project_point_and_pull_back(), can be used as input to this function. If this is not the case, an Exception is thrown.

Definition at line 837 of file utilities.cc.

◆ push_forward_and_differential_forms()

std::tuple< Point< 3 >, Tensor< 1, 3 >, double, double > OpenCASCADE::push_forward_and_differential_forms ( const TopoDS_Face &  face,
const double  u,
const double  v,
const double  tolerance = 1e-7 
)

Given a TopoDS_Face face and the reference coordinates within this face, returns the corresponding point in real space, the normal to the surface at that point and the min and max curvatures as a tuple.

Definition at line 858 of file utilities.cc.

◆ closest_point_and_differential_forms()

std::tuple< Point< 3 >, Tensor< 1, 3 >, double, double > OpenCASCADE::closest_point_and_differential_forms ( const TopoDS_Shape &  in_shape,
const Point< 3 > &  origin,
const double  tolerance = 1e-7 
)

Get the closest point to the given topological shape, together with the normal and the min and max curvatures at that point. If the shape is not elementary, all its sub-faces (only the faces) are iterated, faces first, and only the closest point is returned. This function will throw an exception if the in_shape does not contain at least one face.

Definition at line 801 of file utilities.cc.

◆ line_intersection()

template<int dim>
Point< dim > OpenCASCADE::line_intersection ( const TopoDS_Shape &  in_shape,
const Point< dim > &  origin,
const Tensor< 1, dim > &  direction,
const double  tolerance = 1e-7 
)

Intersect a line passing through the given origin point along direction and the given topological shape. If there is more than one intersection, it will return the closest one.

The optional tolerance parameter is used to compute distances.

Definition at line 511 of file utilities.cc.

◆ point() [1/2]

template<int spacedim>
Point< spacedim > OpenCASCADE::point ( const gp_Pnt &  p,
const double  tolerance = 1e-10 
)

Convert OpenCASCADE point into a Point<spacedim>.

The tolerance argument is used to check if the non used components of the OpenCASCADE point are close to zero. If this is not the case, an assertion is thrown in debug mode.

Definition at line 189 of file utilities.cc.

◆ point() [2/2]

template<int spacedim>
gp_Pnt OpenCASCADE::point ( const Point< spacedim > &  p)

Convert Point<3> into OpenCASCADE point.

Definition at line 172 of file utilities.cc.

◆ point_compare()

template<int dim>
bool OpenCASCADE::point_compare ( const Point< dim > &  p1,
const Point< dim > &  p2,
const Tensor< 1, dim > &  direction = Tensor<1, dim>(),
const double  tolerance = 1e-10 
)

Sort two points according to their scalar product with direction. If the norm of the direction is zero, then use lexicographical ordering. The optional parameter is used as a relative tolerance when comparing objects.

Definition at line 216 of file utilities.cc.