Reference documentation for deal.II version 9.3.3
\(\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\}}\)
Namespaces | Classes | Functions | Variables
OpenCASCADE
Collaboration diagram for OpenCASCADE:

Namespaces

namespace  OpenCASCADE
 

Classes

class  OpenCASCADE::NormalProjectionManifold< dim, spacedim >
 
class  OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >
 
class  OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >
 
class  OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >
 
class  OpenCASCADE::NURBSPatchManifold< dim, spacedim >
 

Functions

 OpenCASCADE::NormalProjectionManifold< dim, spacedim >::NormalProjectionManifold (const TopoDS_Shape &sh, const double tolerance=1e-7)
 
virtual std::unique_ptr< Manifold< dim, spacedim > > OpenCASCADE::NormalProjectionManifold< dim, spacedim >::clone () const override
 
virtual Point< spacedim > OpenCASCADE::NormalProjectionManifold< dim, spacedim >::project_to_manifold (const ArrayView< const Point< spacedim > > &surrounding_points, const Point< spacedim > &candidate) const override
 
 OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::DirectionalProjectionManifold (const TopoDS_Shape &sh, const Tensor< 1, spacedim > &direction, const double tolerance=1e-7)
 
virtual std::unique_ptr< Manifold< dim, spacedim > > OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::clone () const override
 
virtual Point< spacedim > OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::project_to_manifold (const ArrayView< const Point< spacedim > > &surrounding_points, const Point< spacedim > &candidate) const override
 
 OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >::NormalToMeshProjectionManifold (const TopoDS_Shape &sh, const double tolerance=1e-7)
 
virtual std::unique_ptr< Manifold< dim, spacedim > > OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >::clone () const override
 
virtual Point< spacedim > OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >::project_to_manifold (const ArrayView< const Point< spacedim > > &surrounding_points, const Point< spacedim > &candidate) const override
 
 OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::ArclengthProjectionLineManifold (const TopoDS_Shape &sh, const double tolerance=1e-7)
 
virtual std::unique_ptr< Manifold< dim, spacedim > > OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::clone () const override
 
virtual Point< 1 > OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::pull_back (const Point< spacedim > &space_point) const override
 
virtual Point< spacedim > OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::push_forward (const Point< 1 > &chart_point) const override
 
 OpenCASCADE::NURBSPatchManifold< dim, spacedim >::NURBSPatchManifold (const TopoDS_Face &face, const double tolerance=1e-7)
 
virtual std::unique_ptr< Manifold< dim, spacedim > > OpenCASCADE::NURBSPatchManifold< dim, spacedim >::clone () const override
 
virtual Point< 2 > OpenCASCADE::NURBSPatchManifold< dim, spacedim >::pull_back (const Point< spacedim > &space_point) const override
 
virtual Point< spacedim > OpenCASCADE::NURBSPatchManifold< dim, spacedim >::push_forward (const Point< 2 > &chart_point) const override
 
virtual DerivativeForm< 1, 2, spacedim > OpenCASCADE::NURBSPatchManifold< dim, spacedim >::push_forward_gradient (const Point< 2 > &chart_point) const override
 
std::tuple< double, double, double, double > OpenCASCADE::NURBSPatchManifold< dim, spacedim >::get_uv_bounds () const
 

Variables

const TopoDS_Shape OpenCASCADE::NormalProjectionManifold< dim, spacedim >::sh
 
const double OpenCASCADE::NormalProjectionManifold< dim, spacedim >::tolerance
 
const TopoDS_Shape OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::sh
 
const Tensor< 1, spacedim > OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::direction
 
const double OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::tolerance
 
const TopoDS_Shape OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >::sh
 
const double OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >::tolerance
 
const TopoDS_Shape OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::sh
 
Handle_Adaptor3d_HCurve OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::curve
 
const double OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::tolerance
 
const double OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::length
 
TopoDS_Face OpenCASCADE::NURBSPatchManifold< dim, spacedim >::face
 
double OpenCASCADE::NURBSPatchManifold< dim, spacedim >::tolerance
 

Detailed Description

Function Documentation

◆ NormalProjectionManifold()

template<int dim, int spacedim>
OpenCASCADE::NormalProjectionManifold< dim, spacedim >::NormalProjectionManifold ( const TopoDS_Shape &  sh,
const double  tolerance = 1e-7 
)

The standard constructor takes a generic TopoDS_Shape sh, and a tolerance used to compute distances internally.

The TopoDS_Shape can be arbitrary, i.e., a collection of shapes, faces, edges or a single face or edge.

Definition at line 80 of file manifold_lib.cc.

◆ clone() [1/5]

template<int dim, int spacedim>
std::unique_ptr< Manifold< dim, spacedim > > OpenCASCADE::NormalProjectionManifold< dim, spacedim >::clone
overridevirtual

Clone the current Manifold.

Reimplemented from FlatManifold< dim, spacedim >.

Definition at line 93 of file manifold_lib.cc.

◆ project_to_manifold() [1/3]

template<int dim, int spacedim>
Point< spacedim > OpenCASCADE::NormalProjectionManifold< dim, spacedim >::project_to_manifold ( const ArrayView< const Point< spacedim > > &  surrounding_points,
const Point< spacedim > &  candidate 
) const
overridevirtual

Perform the actual projection onto the manifold. This function, in debug mode, checks that each of the surrounding_points is within tolerance from the given TopoDS_Shape. If this is not the case, an exception is thrown.

The projected point is computed using OpenCASCADE normal projection algorithms.

Reimplemented from FlatManifold< dim, spacedim >.

Definition at line 103 of file manifold_lib.cc.

◆ DirectionalProjectionManifold()

template<int dim, int spacedim>
OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::DirectionalProjectionManifold ( const TopoDS_Shape &  sh,
const Tensor< 1, spacedim > &  direction,
const double  tolerance = 1e-7 
)

Construct a Manifold object which will project points on the TopoDS_Shape sh, along the given direction.

Definition at line 121 of file manifold_lib.cc.

◆ clone() [2/5]

template<int dim, int spacedim>
std::unique_ptr< Manifold< dim, spacedim > > OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::clone
overridevirtual

Clone the current Manifold.

Reimplemented from FlatManifold< dim, spacedim >.

Definition at line 136 of file manifold_lib.cc.

◆ project_to_manifold() [2/3]

template<int dim, int spacedim>
Point< spacedim > OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::project_to_manifold ( const ArrayView< const Point< spacedim > > &  surrounding_points,
const Point< spacedim > &  candidate 
) const
overridevirtual

Perform the actual projection onto the manifold. This function, in debug mode, checks that each of the surrounding_points is within tolerance from the given TopoDS_Shape. If this is not the case, an exception is thrown.

The projected point is computed using OpenCASCADE directional projection algorithms.

Reimplemented from FlatManifold< dim, spacedim >.

Definition at line 146 of file manifold_lib.cc.

◆ NormalToMeshProjectionManifold()

template<int dim, int spacedim>
OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >::NormalToMeshProjectionManifold ( const TopoDS_Shape &  sh,
const double  tolerance = 1e-7 
)

Construct a Manifold object which will project points on the TopoDS_Shape sh, along a direction which is approximately normal to the mesh cell.

Definition at line 165 of file manifold_lib.cc.

◆ clone() [3/5]

template<int dim, int spacedim>
std::unique_ptr< Manifold< dim, spacedim > > OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >::clone
overridevirtual

Clone the current Manifold.

Reimplemented from FlatManifold< dim, spacedim >.

Definition at line 180 of file manifold_lib.cc.

◆ project_to_manifold() [3/3]

template<int dim, int spacedim>
Point< spacedim > OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >::project_to_manifold ( const ArrayView< const Point< spacedim > > &  surrounding_points,
const Point< spacedim > &  candidate 
) const
overridevirtual

Perform the actual projection onto the manifold. This function, in debug mode, checks that each of the surrounding_points is within tolerance from the given TopoDS_Shape. If this is not the case, an exception is thrown.

Reimplemented from FlatManifold< dim, spacedim >.

Definition at line 368 of file manifold_lib.cc.

◆ ArclengthProjectionLineManifold()

template<int dim, int spacedim>
OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::ArclengthProjectionLineManifold ( const TopoDS_Shape &  sh,
const double  tolerance = 1e-7 
)

Default constructor with a TopoDS_Edge.

Definition at line 381 of file manifold_lib.cc.

◆ clone() [4/5]

template<int dim, int spacedim>
std::unique_ptr< Manifold< dim, spacedim > > OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::clone
overridevirtual

Clone the current Manifold.

Implements Manifold< dim, spacedim >.

Definition at line 400 of file manifold_lib.cc.

◆ pull_back() [1/2]

template<int dim, int spacedim>
Point< 1 > OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::pull_back ( const Point< spacedim > &  space_point) const
overridevirtual

Given a point on real space, find its arclength parameter. Throws an error in debug mode, if the point is not on the TopoDS_Edge given at construction time.

Implements ChartManifold< dim, spacedim, 1 >.

Definition at line 410 of file manifold_lib.cc.

◆ push_forward() [1/2]

template<int dim, int spacedim>
Point< spacedim > OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::push_forward ( const Point< 1 > &  chart_point) const
overridevirtual

Given an arclength parameter, find its image in real space.

Definition at line 429 of file manifold_lib.cc.

◆ NURBSPatchManifold()

template<int dim, int spacedim>
OpenCASCADE::NURBSPatchManifold< dim, spacedim >::NURBSPatchManifold ( const TopoDS_Face &  face,
const double  tolerance = 1e-7 
)

The constructor takes an OpenCASCADE TopoDS_Face face and an optional tolerance. This class uses the interval OpenCASCADE variables u, v to describe the manifold.

Definition at line 440 of file manifold_lib.cc.

◆ clone() [5/5]

template<int dim, int spacedim>
std::unique_ptr< Manifold< dim, spacedim > > OpenCASCADE::NURBSPatchManifold< dim, spacedim >::clone
overridevirtual

Clone the current Manifold.

Implements Manifold< dim, spacedim >.

Definition at line 450 of file manifold_lib.cc.

◆ pull_back() [2/2]

template<int dim, int spacedim>
Point< 2 > OpenCASCADE::NURBSPatchManifold< dim, spacedim >::pull_back ( const Point< spacedim > &  space_point) const
overridevirtual

Pull back the given point from the Euclidean space. Will return the uv coordinates associated with the point space_point.

Implements ChartManifold< dim, spacedim, 2 >.

Definition at line 460 of file manifold_lib.cc.

◆ push_forward() [2/2]

template<int dim, int spacedim>
Point< spacedim > OpenCASCADE::NURBSPatchManifold< dim, spacedim >::push_forward ( const Point< 2 > &  chart_point) const
overridevirtual

Given a chart_point in the uv coordinate system, this method returns the Euclidean coordinates associated.

Definition at line 476 of file manifold_lib.cc.

◆ push_forward_gradient()

template<int dim, int spacedim>
DerivativeForm< 1, 2, spacedim > OpenCASCADE::NURBSPatchManifold< dim, spacedim >::push_forward_gradient ( const Point< 2 > &  chart_point) const
overridevirtual

Given a point in the spacedim dimensional Euclidean space, this method returns the derivatives of the function \(F\) that maps from the uv coordinate system to the Euclidean coordinate system. In other words, it is a matrix of size \(\text{spacedim}\times\text{chartdim}\).

This function is used in the computations required by the get_tangent_vector() function.

Refer to the general documentation of this class for more information.

Definition at line 486 of file manifold_lib.cc.

◆ get_uv_bounds()

template<int dim, int spacedim>
std::tuple< double, double, double, double > OpenCASCADE::NURBSPatchManifold< dim, spacedim >::get_uv_bounds
protected

Return a tuple representing the minimum and maximum values of u and v. Precisely, it returns (u_min, u_max, v_min, v_max)

Definition at line 517 of file manifold_lib.cc.

Variable Documentation

◆ sh [1/4]

template<int dim, int spacedim>
const TopoDS_Shape OpenCASCADE::NormalProjectionManifold< dim, spacedim >::sh
protected

The topological shape which is used internally to project points. You can construct such a shape by calling the OpenCASCADE::read_IGES() function, which will create a TopoDS_Shape with the geometry contained in the IGES file.

Definition at line 105 of file manifold_lib.h.

◆ tolerance [1/5]

template<int dim, int spacedim>
const double OpenCASCADE::NormalProjectionManifold< dim, spacedim >::tolerance
protected

Relative tolerance used by this class to compute distances.

Definition at line 110 of file manifold_lib.h.

◆ sh [2/4]

template<int dim, int spacedim>
const TopoDS_Shape OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::sh
protected

The topological shape which is used internally to project points. You can construct such a shape by calling the OpenCASCADE::read_IGES() function, which will create a TopoDS_Shape with the geometry contained in the IGES file.

Definition at line 171 of file manifold_lib.h.

◆ direction

template<int dim, int spacedim>
const Tensor<1, spacedim> OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::direction
protected

Direction used to project new points on the shape.

Definition at line 176 of file manifold_lib.h.

◆ tolerance [2/5]

template<int dim, int spacedim>
const double OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >::tolerance
protected

Relative tolerance used by this class to compute distances.

Definition at line 181 of file manifold_lib.h.

◆ sh [3/4]

template<int dim, int spacedim>
const TopoDS_Shape OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >::sh
protected

The topological shape which is used internally to project points. You can construct such a shape by calling the OpenCASCADE::read_IGES() function, which will create a TopoDS_Shape with the geometry contained in the IGES file.

Definition at line 263 of file manifold_lib.h.

◆ tolerance [3/5]

template<int dim, int spacedim>
const double OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >::tolerance
protected

Relative tolerance used by this class to compute distances.

Definition at line 268 of file manifold_lib.h.

◆ sh [4/4]

template<int dim, int spacedim>
const TopoDS_Shape OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::sh
protected

The actual shape used to build this object.

Definition at line 322 of file manifold_lib.h.

◆ curve

template<int dim, int spacedim>
Handle_Adaptor3d_HCurve OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::curve
protected

A Curve adaptor. This is the one which is used in the computations, and it points to the right one above.

Definition at line 328 of file manifold_lib.h.

◆ tolerance [4/5]

template<int dim, int spacedim>
const double OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::tolerance
protected

Relative tolerance used in all internal computations.

Definition at line 333 of file manifold_lib.h.

◆ length

template<int dim, int spacedim>
const double OpenCASCADE::ArclengthProjectionLineManifold< dim, spacedim >::length
protected

The total length of the curve. This is also used as a period if the edge is periodic.

Definition at line 339 of file manifold_lib.h.

◆ face

template<int dim, int spacedim>
TopoDS_Face OpenCASCADE::NURBSPatchManifold< dim, spacedim >::face
protected

An OpenCASCADE TopoDS_Face face given by the CAD.

Definition at line 404 of file manifold_lib.h.

◆ tolerance [5/5]

template<int dim, int spacedim>
double OpenCASCADE::NURBSPatchManifold< dim, spacedim >::tolerance
protected

Tolerance used by OpenCASCADE to identify points in each operation.

Definition at line 410 of file manifold_lib.h.