Reference documentation for deal.II version 9.0.0
|
#include <deal.II/opencascade/boundary_lib.h>
Public Member Functions | |
NormalToMeshProjectionBoundary (const TopoDS_Shape &sh, 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 () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (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 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 ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
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 using OpenCASCADE utilities onto the manifold along a direction which is an estimation of the surrounding points (hence mesh cell) normal.
The direction normal to the mesh is particularly useful because it is the direction in which the mesh is missing nodes. For instance, during the refinement of a cell a new node is initially created around the baricenter of the cell. This location somehow ensures a uniform distance from the nodes of the old cell. Projecting such cell baricenter onto the CAD surface in the direction normal to the original cell will then retain uniform distance from the points of the original cell. Of course, at the stage of mesh generation, no dof handler nor finite element are defined, and such direction has to be estimated. For the case in which 8 surrounding points are present, 4 different triangles are identified with the points assigned, and the normals of such triangles are averaged to obtain the approximation of the normal to the cell.
The case in which 2 surrounding points are present (i.e.:a cell edge is being refined) is of course more tricky. The average of the CAD surface normals at the 2 surrounding points is first computed, and then projected onto the plane normal to the segment linking the surrounding points. This again is an attempt to have the new point with equal distance with respect to the surrounding points
This class only operates with CAD faces and makes the assumption that the shape you pass to it contains at least one face. If that is not the case, an Exception is thrown. 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 normal direction estimated from the surrounding points does not intersect the shape. An exception is thrown when this happens.
Definition at line 230 of file boundary_lib.h.
OpenCASCADE::NormalToMeshProjectionBoundary< dim, spacedim >::NormalToMeshProjectionBoundary | ( | 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 160 of file boundary_lib.cc.
|
overridevirtual |
Clone the current Manifold.
Reimplemented from FlatManifold< dim, spacedim >.
Definition at line 172 of file boundary_lib.cc.
|
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 181 of file boundary_lib.cc.
|
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 263 of file boundary_lib.h.
|
private |
Relative tolerance used by this class to compute distances.
Definition at line 268 of file boundary_lib.h.