Reference documentation for deal.II version 9.0.0
Public Member Functions | Private Attributes | List of all members
OpenCASCADE::DirectionalProjectionBoundary< dim, spacedim > Class Template Reference

#include <deal.II/opencascade/boundary_lib.h>

Inheritance diagram for OpenCASCADE::DirectionalProjectionBoundary< dim, spacedim >:
[legend]

Public Member Functions

 DirectionalProjectionBoundary (const TopoDS_Shape &sh, const Tensor< 1, spacedim > &direction, const double tolerance=1e-7)
 
virtual std::unique_ptr< Manifold< dim, spacedim > > clone () const override
 
virtual Point< spacedim > project_to_manifold (const ArrayView< const Point< spacedim >> &surrounding_points, const Point< spacedim > &candidate) const override
 
- Public Member Functions inherited from FlatManifold< dim, spacedim >
 FlatManifold (const Tensor< 1, spacedim > &periodicity=Tensor< 1, spacedim >(), const double tolerance=1e-10)
 
virtual Point< spacedim > get_new_point (const ArrayView< const Point< spacedim >> &surrounding_points, const ArrayView< const double > &weights) const override
 
virtual void get_new_points (const ArrayView< const Point< spacedim >> &surrounding_points, const Table< 2, double > &weights, ArrayView< Point< spacedim >> new_points) const override
 
virtual Tensor< 1, spacedim > get_tangent_vector (const Point< spacedim > &x1, const Point< spacedim > &x2) const override
 
virtual Tensor< 1, spacedim > normal_vector (const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const override
 
virtual void get_normals_at_vertices (const typename Triangulation< dim, spacedim >::face_iterator &face, typename Manifold< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const override
 
const Tensor< 1, spacedim > & get_periodicity () const
 
- Public Member Functions inherited from Manifold< dim, spacedim >
virtual ~Manifold ()=default
 
virtual Point< spacedim > get_intermediate_point (const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const
 
virtual Point< spacedim > get_new_point_on_line (const typename Triangulation< dim, spacedim >::line_iterator &line) const
 
virtual Point< spacedim > get_new_point_on_quad (const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
 
virtual Point< spacedim > get_new_point_on_hex (const typename Triangulation< dim, spacedim >::hex_iterator &hex) const
 
Point< spacedim > get_new_point_on_face (const typename Triangulation< dim, spacedim >::face_iterator &face) const
 
Point< spacedim > get_new_point_on_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
 
virtual void get_normals_at_vertices (const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (const char *identifier=nullptr) const
 
void unsubscribe (const char *identifier=nullptr) const
 
unsigned int n_subscriptions () const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Private Attributes

const TopoDS_Shape sh
 
const Tensor< 1, spacedim > direction
 
const double tolerance
 

Additional Inherited Members

- Public Types inherited from Manifold< dim, spacedim >
typedef Tensor< 1, spacedim > FaceVertexNormals[GeometryInfo< dim >::vertices_per_face]
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Detailed Description

template<int dim, int spacedim>
class OpenCASCADE::DirectionalProjectionBoundary< dim, spacedim >

A Manifold object based on OpenCASCADE TopoDS_Shape where new points are first computed by averaging the surrounding points in the same way as FlatManifold does, and then projecting them onto the manifold along the direction specified at construction time using OpenCASCADE utilities.

This class makes no assumptions on the shape you pass to it, and the topological dimension of the Manifold is inferred from the TopoDS_Shape itself. In debug mode there is a sanity check to make sure that the surrounding points (the ones used in project_to_manifold()) actually live on the Manifold, i.e., calling OpenCASCADE::closest_point() on those points leaves them untouched. If this is not the case, an ExcPointNotOnManifold is thrown.

Notice that this type of Manifold descriptor may fail to give results if the triangulation to be refined is close to the boundary of the given TopoDS_Shape, or when the direction you use at construction time does not intersect the shape. An exception is thrown when this happens.

Author
Luca Heltai, Andrea Mola, 2011–2014.

Definition at line 135 of file boundary_lib.h.

Constructor & Destructor Documentation

◆ DirectionalProjectionBoundary()

template<int dim, int spacedim>
OpenCASCADE::DirectionalProjectionBoundary< dim, spacedim >::DirectionalProjectionBoundary ( 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 118 of file boundary_lib.cc.

Member Function Documentation

◆ clone()

template<int dim, int spacedim>
std::unique_ptr< Manifold< dim, spacedim > > OpenCASCADE::DirectionalProjectionBoundary< dim, spacedim >::clone ( ) const
overridevirtual

Clone the current Manifold.

Reimplemented from FlatManifold< dim, spacedim >.

Definition at line 132 of file boundary_lib.cc.

◆ project_to_manifold()

template<int dim, int spacedim>
Point< spacedim > OpenCASCADE::DirectionalProjectionBoundary< 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 142 of file boundary_lib.cc.

Member Data Documentation

◆ sh

template<int dim, int spacedim>
const TopoDS_Shape OpenCASCADE::DirectionalProjectionBoundary< dim, spacedim >::sh
private

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 boundary_lib.h.

◆ direction

template<int dim, int spacedim>
const Tensor<1, spacedim> OpenCASCADE::DirectionalProjectionBoundary< dim, spacedim >::direction
private

Direction used to project new points on the shape.

Definition at line 176 of file boundary_lib.h.

◆ tolerance

template<int dim, int spacedim>
const double OpenCASCADE::DirectionalProjectionBoundary< dim, spacedim >::tolerance
private

Relative tolerance used by this class to compute distances.

Definition at line 181 of file boundary_lib.h.


The documentation for this class was generated from the following files: