dolfin/geometry¶
Documentation for C++ code found in dolfin/geometry/*.h
Contents
Functions¶
intersect¶
C++ documentation for intersect
from dolfin/geometry/intersect.h
:
-
std::shared_ptr<const MeshPointIntersection>
dolfin::
intersect
(const Mesh &mesh, const Point &point)¶ Compute and return intersection between
Mesh
andPoint
. Arguments mesh (Mesh
) The mesh to be intersected. point (Point
) The point to be intersected. Returns:cpp:any:MeshPointIntersection The intersection data.Parameters: - mesh –
- point –
operator<<¶
C++ documentation for operator<<
from dolfin/geometry/Point.h
:
Classes¶
BoundingBoxTree¶
C++ documentation for BoundingBoxTree
from dolfin/geometry/BoundingBoxTree.h
:
-
class
dolfin::
BoundingBoxTree
¶ This class implements a (distributed) axis aligned bounding box tree (AABB tree). Bounding box trees can be created from meshes and [other data structures, to be filled in].
-
dolfin::BoundingBoxTree::
BoundingBoxTree
()¶ Create empty bounding box tree.
-
void
dolfin::BoundingBoxTree::
build
(const Mesh &mesh)¶ Build bounding box tree for cells of mesh. Arguments mesh (
Mesh
) The mesh for which to compute the bounding box tree.Parameters: mesh –
-
void
dolfin::BoundingBoxTree::
build
(const Mesh &mesh, std::size_t tdim)¶ Build bounding box tree for mesh entities of given dimension. Arguments mesh (
Mesh
) The mesh for which to compute the bounding box tree. dimension (std::size_t) The entity dimension (topological dimension) for which to compute the bounding box tree.Parameters: - mesh –
- tdim –
-
void
dolfin::BoundingBoxTree::
build
(const std::vector<Point> &points, std::size_t gdim)¶ Build bounding box tree for point cloud. Arguments points (std::vector<
Point
>) The list of points. gdim (std::size_t) The geometric dimension.Parameters: - points –
- gdim –
-
bool
dolfin::BoundingBoxTree::
collides
(const Point &point) const¶ Check whether given point collides with the bounding box tree. This is equivalent to calling compute_first_collision and checking whether any collision was detected. Returns bool True iff the point is inside the tree.
Parameters: point –
-
bool
dolfin::BoundingBoxTree::
collides_entity
(const Point &point) const¶ Check whether given point collides with any entity contained in the bounding box tree. This is equivalent to calling compute_first_entity_collision and checking whether any collision was detected. Returns bool True iff the point is inside the tree.
Parameters: point –
-
std::pair<unsigned int, double>
dolfin::BoundingBoxTree::
compute_closest_entity
(const Point &point) const¶ Compute closest entity to
Point
. Returns unsigned int The local index for the entity that is closest to the point. If more than one entity is at the same distance (or point contained in entity), then the first entity is returned. double The distance to the closest entity. Arguments point (Point
) The point.Parameters: point –
-
std::pair<unsigned int, double>
dolfin::BoundingBoxTree::
compute_closest_point
(const Point &point) const¶ Compute closest point to
Point
. This function assumes that the tree has been built for a point cloud. Developer note: This function should not be confused with computing the closest point in all entities of a mesh. That function could be added with relative ease since we actually compute the closest points to get the distance in the above function (compute_closest_entity) inside the specialized implementations in TetrahedronCell.cpp etc. Returns unsigned int The local index for the point that is closest to the point. If more than one point is at the same distance (or point contained in entity), then the first point is returned. double The distance to the closest point. Arguments point (Point
) The point.Parameters: point –
-
std::pair<std::vector<unsigned int>, std::vector<unsigned int>>
dolfin::BoundingBoxTree::
compute_collisions
(const BoundingBoxTree &tree) const¶ Compute all collisions between bounding boxes and
BoundingBoxTree
. Returns std::vector<unsigned int> A list of local indices for entities in this tree that collide with (intersect) entities in other tree. std::vector<unsigned int> A list of local indices for entities in other tree that collide with (intersect) entities in this tree. The two lists have equal length and contain matching entities, such that entityi
in the first list collides with entityi
in the second list. Note that this means that the entity lists may contain duplicate entities since a single entity may collide with several different entities. Arguments tree (BoundingBoxTree
) The bounding box tree. Note that this function only checks collisions between bounding boxes of entities. It does not check that the entities themselves actually collide. To compute entity collisions, use the function compute_entity_collisions.Parameters: tree –
-
std::vector<unsigned int>
dolfin::BoundingBoxTree::
compute_collisions
(const Point &point) const¶ Compute all collisions between bounding boxes and
Point
. Returns std::vector<unsigned int> A list of local indices for entities contained in (leaf) bounding boxes that collide with (intersect) the given point. Arguments point (Point
) The point.Parameters: point –
-
std::pair<std::vector<unsigned int>, std::vector<unsigned int>>
dolfin::BoundingBoxTree::
compute_entity_collisions
(const BoundingBoxTree &tree) const¶ Compute all collisions between entities and
BoundingBoxTree
. Returns std::vector<unsigned int> A list of local indices for entities in this tree that collide with (intersect) entities in other tree. std::vector<unsigned int> A list of local indices for entities in other tree that collide with (intersect) entities in this tree. The two lists have equal length and contain matching entities, such that entityi
in the first list collides with entityi
in the second list. Note that this means that the entity lists may contain duplicate entities since a single entity may collide with several different entities. Arguments tree (BoundingBoxTree
) The bounding box tree.Parameters: tree –
-
std::vector<unsigned int>
dolfin::BoundingBoxTree::
compute_entity_collisions
(const Point &point) const¶ Compute all collisions between entities and
Point
. Returns std::vector<unsigned int> A list of local indices for entities that collide with (intersect) the given point. Arguments point (Point
) The point.Parameters: point –
-
unsigned int
dolfin::BoundingBoxTree::
compute_first_collision
(const Point &point) const¶ Compute first collision between bounding boxes and
Point
. Returns unsigned int The local index for the first found entity contained in a (leaf) bounding box that collides with (intersects) the given point. If not found, std::numeric_limits<unsigned int>::max() is returned. Arguments point (Point
) The point.Parameters: point –
-
unsigned int
dolfin::BoundingBoxTree::
compute_first_entity_collision
(const Point &point) const¶ Compute first collision between entities and
Point
. Returns unsigned int The local index for the first found entity that collides with (intersects) the given point. If not found, std::numeric_limits<unsigned int>::max() is returned. Arguments point (Point
) The point.Parameters: point –
-
std::vector<unsigned int>
dolfin::BoundingBoxTree::
compute_process_collisions
(const Point &point) const¶ Compute all collisions between process bounding boxes and
Point
. Effectively a list of processes which may contain thePoint
. Returns std::vector<unsigned int> A list of process numbers where theMesh
may collide with (intersect) the given point. Arguments point (Point
) The point.Parameters: point –
-
dolfin::BoundingBoxTree::
~BoundingBoxTree
()¶ Destructor.
-
BoundingBoxTree1D¶
C++ documentation for BoundingBoxTree1D
from dolfin/geometry/BoundingBoxTree1D.h
:
-
class
dolfin::
BoundingBoxTree1D
¶ Specialization of bounding box implementation to 1D.
-
bool
dolfin::BoundingBoxTree1D::
bbox_in_bbox
(const double *a, unsigned int node) const¶ Check whether bounding box (a) collides with bounding box (node)
Parameters: - a –
- node –
-
void
dolfin::BoundingBoxTree1D::
compute_bbox_of_bboxes
(double *bbox, std::size_t &axis, const std::vector<double> &leaf_bboxes, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &end)¶ Compute bounding box of bounding boxes.
Parameters: - bbox –
- axis –
- leaf_bboxes –
- begin –
- end –
-
void
dolfin::BoundingBoxTree1D::
compute_bbox_of_points
(double *bbox, std::size_t &axis, const std::vector<Point> &points, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &end)¶ Compute bounding box of points.
Parameters: - bbox –
- axis –
- points –
- begin –
- end –
-
double
dolfin::BoundingBoxTree1D::
compute_squared_distance_bbox
(const double *x, unsigned int node) const¶ Compute squared distance between point and bounding box.
Parameters: - x –
- node –
-
double
dolfin::BoundingBoxTree1D::
compute_squared_distance_point
(const double *x, unsigned int node) const¶ Compute squared distance between point and point.
Parameters: - x –
- node –
-
std::size_t
dolfin::BoundingBoxTree1D::
gdim
() const¶ Return geometric dimension.
-
const double *
dolfin::BoundingBoxTree1D::
get_bbox_coordinates
(unsigned int node) const¶ Return bounding box coordinates for node.
Parameters: node –
-
bool
dolfin::BoundingBoxTree1D::
point_in_bbox
(const double *x, unsigned int node) const¶ Check whether point (x) is in bounding box (node)
Parameters: - x –
- node –
-
void
dolfin::BoundingBoxTree1D::
sort_bboxes
(std::size_t axis, const std::vector<double> &leaf_bboxes, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &middle, const std::vector<unsigned int>::iterator &end)¶ Sort leaf bounding boxes along given axis.
Parameters: - axis –
- leaf_bboxes –
- begin –
- middle –
- end –
-
bool
BoundingBoxTree2D¶
C++ documentation for BoundingBoxTree2D
from dolfin/geometry/BoundingBoxTree2D.h
:
-
class
dolfin::
BoundingBoxTree2D
¶ Specialization of bounding box implementation to 2D.
-
bool
dolfin::BoundingBoxTree2D::
bbox_in_bbox
(const double *a, unsigned int node) const¶ Check whether bounding box (a) collides with bounding box (node)
Parameters: - a –
- node –
-
void
dolfin::BoundingBoxTree2D::
compute_bbox_of_bboxes
(double *bbox, std::size_t &axis, const std::vector<double> &leaf_bboxes, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &end)¶ Compute bounding box of bounding boxes.
Parameters: - bbox –
- axis –
- leaf_bboxes –
- begin –
- end –
-
void
dolfin::BoundingBoxTree2D::
compute_bbox_of_points
(double *bbox, std::size_t &axis, const std::vector<Point> &points, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &end)¶ Compute bounding box of points.
Parameters: - bbox –
- axis –
- points –
- begin –
- end –
-
double
dolfin::BoundingBoxTree2D::
compute_squared_distance_bbox
(const double *x, unsigned int node) const¶ Compute squared distance between point and bounding box.
Parameters: - x –
- node –
-
double
dolfin::BoundingBoxTree2D::
compute_squared_distance_point
(const double *x, unsigned int node) const¶ Compute squared distance between point and point.
Parameters: - x –
- node –
-
std::size_t
dolfin::BoundingBoxTree2D::
gdim
() const¶ Return geometric dimension.
-
const double *
dolfin::BoundingBoxTree2D::
get_bbox_coordinates
(unsigned int node) const¶ Return bounding box coordinates for node.
Parameters: node –
-
bool
dolfin::BoundingBoxTree2D::
point_in_bbox
(const double *x, unsigned int node) const¶ Check whether point (x) is in bounding box (node)
Parameters: - x –
- node –
-
void
dolfin::BoundingBoxTree2D::
sort_bboxes
(std::size_t axis, const std::vector<double> &leaf_bboxes, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &middle, const std::vector<unsigned int>::iterator &end)¶ Sort leaf bounding boxes along given axis.
Parameters: - axis –
- leaf_bboxes –
- begin –
- middle –
- end –
-
bool
BoundingBoxTree3D¶
C++ documentation for BoundingBoxTree3D
from dolfin/geometry/BoundingBoxTree3D.h
:
-
class
dolfin::
BoundingBoxTree3D
¶ Specialization of bounding box implementation to 3D.
-
bool
dolfin::BoundingBoxTree3D::
bbox_in_bbox
(const double *a, unsigned int node) const¶ Check whether bounding box (a) collides with bounding box (node)
Parameters: - a –
- node –
-
void
dolfin::BoundingBoxTree3D::
compute_bbox_of_bboxes
(double *bbox, std::size_t &axis, const std::vector<double> &leaf_bboxes, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &end)¶ Compute bounding box of bounding boxes.
Parameters: - bbox –
- axis –
- leaf_bboxes –
- begin –
- end –
-
void
dolfin::BoundingBoxTree3D::
compute_bbox_of_points
(double *bbox, std::size_t &axis, const std::vector<Point> &points, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &end)¶ Compute bounding box of points.
Parameters: - bbox –
- axis –
- points –
- begin –
- end –
-
double
dolfin::BoundingBoxTree3D::
compute_squared_distance_bbox
(const double *x, unsigned int node) const¶ Compute squared distance between point and bounding box.
Parameters: - x –
- node –
-
double
dolfin::BoundingBoxTree3D::
compute_squared_distance_point
(const double *x, unsigned int node) const¶ Compute squared distance between point and point.
Parameters: - x –
- node –
-
std::size_t
dolfin::BoundingBoxTree3D::
gdim
() const¶ Return geometric dimension.
-
const double *
dolfin::BoundingBoxTree3D::
get_bbox_coordinates
(unsigned int node) const¶ Return bounding box coordinates for node.
Parameters: node –
-
bool
dolfin::BoundingBoxTree3D::
point_in_bbox
(const double *x, const unsigned int node) const¶ Check whether point (x) is in bounding box (node)
Parameters: - x –
- node –
-
void
dolfin::BoundingBoxTree3D::
sort_bboxes
(std::size_t axis, const std::vector<double> &leaf_bboxes, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &middle, const std::vector<unsigned int>::iterator &end)¶ Sort leaf bounding boxes along given axis.
Parameters: - axis –
- leaf_bboxes –
- begin –
- middle –
- end –
-
bool
CollisionDetection¶
C++ documentation for CollisionDetection
from dolfin/geometry/CollisionDetection.h
:
-
class
dolfin::
CollisionDetection
¶ This class implements algorithms for detecting pairwise collisions between mesh entities of varying dimensions.
-
bool
dolfin::CollisionDetection::
collides
(const MeshEntity &entity, const Point &point)¶ Check whether entity collides with point.
Parameters: - entity – (
MeshEntity
) The entity. - point – (
Point
) The point.
Returns: bool True iff entity collides with cell.
- entity – (
-
bool
dolfin::CollisionDetection::
collides
(const MeshEntity &entity_0, const MeshEntity &entity_1)¶ Check whether two entities collide.
Parameters: - entity_0 – (
MeshEntity
) The first entity. - entity_1 – (
MeshEntity
) The second entity.
Returns: bool True iff entity collides with cell.
- entity_0 – (
-
bool
dolfin::CollisionDetection::
collides_edge_edge
(const Point &a, const Point &b, const Point &c, const Point &d)¶ Check whether edge a-b collides with edge c-d.
Parameters: - a –
- b –
- c –
- d –
-
bool
dolfin::CollisionDetection::
collides_interval_interval
(const MeshEntity &interval_0, const MeshEntity &interval_1)¶ Check whether interval collides with interval.
Parameters: - interval_0 – (
MeshEntity
) The first interval. - interval_1 – (
MeshEntity
) The second interval.
Returns: bool True iff objects collide.
- interval_0 – (
-
bool
dolfin::CollisionDetection::
collides_interval_point
(const MeshEntity &interval, const Point &point)¶ Check whether interval collides with point.
Parameters: - interval – (
MeshEntity
) The interval. - point – (
Point
) The point.
Returns: bool True iff objects collide.
- interval – (
-
bool
dolfin::CollisionDetection::
collides_interval_point
(const Point &p0, const Point &p1, const Point &point)¶ The implementation of collides_interval_point.
Parameters: - p0 –
- p1 –
- point –
-
bool
dolfin::CollisionDetection::
collides_tetrahedron_point
(const MeshEntity &tetrahedron, const Point &point)¶ Check whether tetrahedron collides with point.
Parameters: - tetrahedron – (
MeshEntity
) The tetrahedron. - point – (
Point
) The point.
Returns: bool True iff objects collide.
- tetrahedron – (
-
bool
dolfin::CollisionDetection::
collides_tetrahedron_point
(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &point)¶ The implementation of collides_tetrahedron_point.
Parameters: - p0 –
- p1 –
- p2 –
- p3 –
- point –
-
bool
dolfin::CollisionDetection::
collides_tetrahedron_tetrahedron
(const MeshEntity &tetrahedron_0, const MeshEntity &tetrahedron_1)¶ Check whether tetrahedron collides with tetrahedron.
Parameters: - tetrahedron_0 – (
MeshEntity
) The first tetrahedron. - tetrahedron_1 – (
MeshEntity
) The second tetrahedron.
Returns: bool True iff objects collide.
- tetrahedron_0 – (
-
bool
dolfin::CollisionDetection::
collides_tetrahedron_triangle
(const MeshEntity &tetrahedron, const MeshEntity &triangle)¶ Check whether tetrahedron collides with triangle.
Parameters: - tetrahedron – (
MeshEntity
) The tetrahedron. - triangle – (
MeshEntity
) The triangle.
Returns: bool True iff objects collide.
- tetrahedron – (
-
bool
dolfin::CollisionDetection::
collides_tetrahedron_triangle
(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2)¶ Parameters: - p0 –
- p1 –
- p2 –
- p3 –
- q0 –
- q1 –
- q2 –
-
bool
dolfin::CollisionDetection::
collides_triangle_point
(const MeshEntity &triangle, const Point &point)¶ Check whether triangle collides with point.
Parameters: - triangle – (
MeshEntity
) The triangle. - point – (
Point
) The point.
Returns: bool True iff objects collide.
- triangle – (
-
bool
dolfin::CollisionDetection::
collides_triangle_point
(const Point &p0, const Point &p1, const Point &p2, const Point &point)¶ The implementation of collides_triangle_point.
Parameters: - p0 –
- p1 –
- p2 –
- point –
-
bool
dolfin::CollisionDetection::
collides_triangle_point_2d
(const Point &p0, const Point &p1, const Point &p2, const Point &point)¶ Specialised implementation of collides_triangle_point in 2D.
Parameters: - p0 –
- p1 –
- p2 –
- point –
-
bool
dolfin::CollisionDetection::
collides_triangle_triangle
(const MeshEntity &triangle_0, const MeshEntity &triangle_1)¶ Check whether triangle collides with triangle.
Parameters: - triangle_0 – (
MeshEntity
) The first triangle. - triangle_1 – (
MeshEntity
) The second triangle.
Returns: bool True iff objects collide.
- triangle_0 – (
-
bool
dolfin::CollisionDetection::
collides_triangle_triangle
(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)¶ Parameters: - p0 –
- p1 –
- p2 –
- q0 –
- q1 –
- q2 –
-
bool
dolfin::CollisionDetection::
compute_intervals
(double VV0, double VV1, double VV2, double D0, double D1, double D2, double D0D1, double D0D2, double &A, double &B, double &C, double &X0, double &X1)¶ Parameters: - VV0 –
- VV1 –
- VV2 –
- D0 –
- D1 –
- D2 –
- D0D1 –
- D0D2 –
- A –
- B –
- C –
- X0 –
- X1 –
-
bool
dolfin::CollisionDetection::
coplanar_tri_tri
(const Point &N, const Point &V0, const Point &V1, const Point &V2, const Point &U0, const Point &U1, const Point &U2)¶ Parameters: - N –
- V0 –
- V1 –
- V2 –
- U0 –
- U1 –
- U2 –
-
bool
dolfin::CollisionDetection::
edge_against_tri_edges
(int i0, int i1, const Point &V0, const Point &V1, const Point &U0, const Point &U1, const Point &U2)¶ Parameters: - i0 –
- i1 –
- V0 –
- V1 –
- U0 –
- U1 –
- U2 –
-
bool
dolfin::CollisionDetection::
edge_edge_test
(int i0, int i1, double Ax, double Ay, const Point &V0, const Point &U0, const Point &U1)¶ Parameters: - i0 –
- i1 –
- Ax –
- Ay –
- V0 –
- U0 –
- U1 –
-
bool
dolfin::CollisionDetection::
point_in_tri
(int i0, int i1, const Point &V0, const Point &U0, const Point &U1, const Point &U2)¶ Parameters: - i0 –
- i1 –
- V0 –
- U0 –
- U1 –
- U2 –
-
bool
dolfin::CollisionDetection::
separating_plane_edge_A
(const std::vector<std::vector<double>> &coord_1, const std::vector<int> &masks, int f0, int f1)¶ Parameters: - coord_1 –
- masks –
- f0 –
- f1 –
-
bool
dolfin::CollisionDetection::
separating_plane_face_A_1
(const std::vector<Point> &pv1, const Point &n, std::vector<double> &bc, int &mask_edges)¶ Parameters: - pv1 –
- n –
- bc –
- mask_edges –
-
bool
dolfin::CollisionDetection::
separating_plane_face_A_2
(const std::vector<Point> &v1, const std::vector<Point> &v2, const Point &n, std::vector<double> &bc, int &mask_edges)¶ Parameters: - v1 –
- v2 –
- n –
- bc –
- mask_edges –
-
bool
GenericBoundingBoxTree¶
C++ documentation for GenericBoundingBoxTree
from dolfin/geometry/GenericBoundingBoxTree.h
:
-
class
dolfin::
GenericBoundingBoxTree
¶ Base class for bounding box implementations (envelope-letter design)
-
dolfin::GenericBoundingBoxTree::
GenericBoundingBoxTree
()¶ Constructor.
-
unsigned int
dolfin::GenericBoundingBoxTree::
add_bbox
(const BBox &bbox, const double *b, std::size_t gdim)¶ Add bounding box and coordinates.
Parameters: - bbox –
- b –
- gdim –
-
unsigned int
dolfin::GenericBoundingBoxTree::
add_point
(const BBox &bbox, const Point &point, std::size_t gdim)¶ Add bounding box and point coordinates.
Parameters: - bbox –
- point –
- gdim –
-
bool
dolfin::GenericBoundingBoxTree::
bbox_in_bbox
(const double *a, unsigned int node) const = 0¶ Check whether bounding box (a) collides with bounding box (node)
Parameters: - a –
- node –
-
void
dolfin::GenericBoundingBoxTree::
build
(const Mesh &mesh, std::size_t tdim)¶ Build bounding box tree for mesh entities of given dimension.
Parameters: - mesh –
- tdim –
-
void
dolfin::GenericBoundingBoxTree::
build
(const std::vector<Point> &points)¶ Build bounding box tree for point cloud.
Parameters: points –
-
void
dolfin::GenericBoundingBoxTree::
build_point_search_tree
(const Mesh &mesh) const¶ Compute point search tree if not already done.
Parameters: mesh –
-
void
dolfin::GenericBoundingBoxTree::
clear
()¶ Clear existing data if any.
-
void
dolfin::GenericBoundingBoxTree::
compute_bbox_of_bboxes
(double *bbox, std::size_t &axis, const std::vector<double> &leaf_bboxes, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &end) = 0¶ Compute bounding box of bounding boxes.
Parameters: - bbox –
- axis –
- leaf_bboxes –
- begin –
- end –
-
void
dolfin::GenericBoundingBoxTree::
compute_bbox_of_entity
(double *b, const MeshEntity &entity, std::size_t gdim) const¶ Compute bounding box of mesh entity.
Parameters: - b –
- entity –
- gdim –
-
void
dolfin::GenericBoundingBoxTree::
compute_bbox_of_points
(double *bbox, std::size_t &axis, const std::vector<Point> &points, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &end) = 0¶ Compute bounding box of points.
Parameters: - bbox –
- axis –
- points –
- begin –
- end –
-
std::pair<unsigned int, double>
dolfin::GenericBoundingBoxTree::
compute_closest_entity
(const Point &point, const Mesh &mesh) const¶ Compute closest entity and distance to
Point
Parameters: - point –
- mesh –
-
std::pair<unsigned int, double>
dolfin::GenericBoundingBoxTree::
compute_closest_point
(const Point &point) const¶ Compute closest point and distance to
Point
Parameters: point –
-
std::pair<std::vector<unsigned int>, std::vector<unsigned int>>
dolfin::GenericBoundingBoxTree::
compute_collisions
(const GenericBoundingBoxTree &tree) const¶ Compute all collisions between bounding boxes and
BoundingBoxTree
Parameters: tree –
-
std::vector<unsigned int>
dolfin::GenericBoundingBoxTree::
compute_collisions
(const Point &point) const¶ Compute all collisions between bounding boxes and
Point
Parameters: point –
-
std::pair<std::vector<unsigned int>, std::vector<unsigned int>>
dolfin::GenericBoundingBoxTree::
compute_entity_collisions
(const GenericBoundingBoxTree &tree, const Mesh &mesh_A, const Mesh &mesh_B) const¶ Compute all collisions between entities and
BoundingBoxTree
Parameters: - tree –
- mesh_A –
- mesh_B –
-
std::vector<unsigned int>
dolfin::GenericBoundingBoxTree::
compute_entity_collisions
(const Point &point, const Mesh &mesh) const¶ Compute all collisions between entities and
Point
Parameters: - point –
- mesh –
-
unsigned int
dolfin::GenericBoundingBoxTree::
compute_first_collision
(const Point &point) const¶ Compute first collision between bounding boxes and
Point
Parameters: point –
-
unsigned int
dolfin::GenericBoundingBoxTree::
compute_first_entity_collision
(const Point &point, const Mesh &mesh) const¶ Compute first collision between entities and
Point
Parameters: - point –
- mesh –
-
std::vector<unsigned int>
dolfin::GenericBoundingBoxTree::
compute_process_collisions
(const Point &point) const¶ Compute all collisions between processes and
Point
Parameters: point –
-
double
dolfin::GenericBoundingBoxTree::
compute_squared_distance_bbox
(const double *x, unsigned int node) const = 0¶ Compute squared distance between point and bounding box.
Parameters: - x –
- node –
-
double
dolfin::GenericBoundingBoxTree::
compute_squared_distance_point
(const double *x, unsigned int node) const = 0¶ Compute squared distance between point and point.
Parameters: - x –
- node –
-
std::shared_ptr<GenericBoundingBoxTree>
dolfin::GenericBoundingBoxTree::
create
(unsigned int dim)¶ Factory function returning (empty) tree of appropriate dimension.
Parameters: dim –
-
std::size_t
dolfin::GenericBoundingBoxTree::
gdim
() const = 0¶ Return geometric dimension.
-
const BBox &
dolfin::GenericBoundingBoxTree::
get_bbox
(unsigned int node) const¶ Return bounding box for given node.
Parameters: node –
-
const double *
dolfin::GenericBoundingBoxTree::
get_bbox_coordinates
(unsigned int node) const = 0¶ Return bounding box coordinates for node.
Parameters: node –
-
bool
dolfin::GenericBoundingBoxTree::
is_leaf
(const BBox &bbox, unsigned int node) const¶ Check whether bounding box is a leaf node.
Parameters: - bbox –
- node –
-
unsigned int
dolfin::GenericBoundingBoxTree::
num_bboxes
() const¶ Return number of bounding boxes.
-
bool
dolfin::GenericBoundingBoxTree::
point_in_bbox
(const double *x, unsigned int node) const = 0¶ Check whether point (x) is in bounding box (node)
Parameters: - x –
- node –
-
void
dolfin::GenericBoundingBoxTree::
sort_bboxes
(std::size_t axis, const std::vector<double> &leaf_bboxes, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &middle, const std::vector<unsigned int>::iterator &end) = 0¶ Sort leaf bounding boxes along given axis.
Parameters: - axis –
- leaf_bboxes –
- begin –
- middle –
- end –
-
void
dolfin::GenericBoundingBoxTree::
sort_points
(std::size_t axis, const std::vector<Point> &points, const std::vector<unsigned int>::iterator &begin, const std::vector<unsigned int>::iterator &middle, const std::vector<unsigned int>::iterator &end)¶ Sort points along given axis.
Parameters: - axis –
- points –
- begin –
- middle –
- end –
-
std::string
dolfin::GenericBoundingBoxTree::
str
(bool verbose = false)¶ Print out for debugging.
Parameters: verbose –
-
void
dolfin::GenericBoundingBoxTree::
tree_print
(std::stringstream &s, unsigned int i)¶ Print out recursively, for debugging.
Parameters: - s –
- i –
-
dolfin::GenericBoundingBoxTree::
~GenericBoundingBoxTree
()¶ Destructor.
-
IntersectionTriangulation¶
C++ documentation for IntersectionTriangulation
from dolfin/geometry/IntersectionTriangulation.h
:
-
class
dolfin::
IntersectionTriangulation
¶ This class implements algorithms for computing triangulations of pairwise intersections of simplices.
-
bool
dolfin::IntersectionTriangulation::
intersection_edge_edge
(const Point &a, const Point &b, const Point &c, const Point &d, Point &pt)¶ Parameters: - a –
- b –
- c –
- d –
- pt –
-
bool
dolfin::IntersectionTriangulation::
intersection_face_edge
(const Point &r, const Point &s, const Point &t, const Point &a, const Point &b, Point &pt)¶ Parameters: - r –
- s –
- t –
- a –
- b –
- pt –
-
void
dolfin::IntersectionTriangulation::
triangulate_intersection
(const MeshEntity &cell, const std::vector<double> &triangulation, const std::vector<Point> &normals, std::vector<double> &intersection_triangulation, std::vector<Point> &intersection_normals, std::size_t tdim)¶ Function
for computing the intersection of a cell with a flat vector of simplices with topological dimension tdim. The geometrical dimension is assumed to be the same as for the cell. The corresponding normals are also saved.Parameters: - cell –
- triangulation –
- normals –
- intersection_triangulation –
- intersection_normals –
- tdim –
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection
(const MeshEntity &cell, const std::vector<double> &triangulation, std::size_t tdim)¶ Function
for computing the intersection of a cell with a flat vector of simplices with topological dimension tdim. The geometrical dimension is assumed to be the same as for the cell.Parameters: - cell –
- triangulation –
- tdim –
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection
(const MeshEntity &entity_0, const MeshEntity &entity_1)¶ Compute triangulation of intersection of two entities
Parameters: - entity_0 – (
MeshEntity
) The first entity. - entity_1 – (
MeshEntity
) The second entity.
Returns: std::vector<double> A flattened array of simplices of dimension num_simplices x num_vertices x gdim = num_simplices x (tdim + 1) x gdim
- entity_0 – (
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection
(const std::vector<Point> &s0, std::size_t tdim0, const std::vector<Point> &s1, std::size_t tdim1, std::size_t gdim)¶ Function
for general intersection computation of two simplices with different topological dimension but the same geometrical dimensionParameters: - s0 –
- tdim0 –
- s1 –
- tdim1 –
- gdim –
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection_interval_interval
(const MeshEntity &interval_0, const MeshEntity &interval_1)¶ Compute triangulation of intersection of two intervals
Parameters: - interval_0 – (
MeshEntity
) The first interval. - interval_1 – (
MeshEntity
) The second interval.
Returns: std::vector<double> A flattened array of simplices of dimension num_simplices x num_vertices x gdim = num_simplices x (tdim + 1) x gdim
- interval_0 – (
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection_interval_interval
(const std::vector<Point> &interval_0, const std::vector<Point> &interval_1, std::size_t gdim)¶ Parameters: - interval_0 –
- interval_1 –
- gdim –
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection_tetrahedron_tetrahedron
(const MeshEntity &tetrahedron_0, const MeshEntity &tetrahedron_1)¶ Compute triangulation of intersection of two tetrahedra
Parameters: - tetrahedron_0 – (
MeshEntity
) The first tetrahedron. - tetrahedron_1 – (
MeshEntity
) The second tetrahedron.
Returns: std::vector<double> A flattened array of simplices of dimension num_simplices x num_vertices x gdim = num_simplices x (tdim + 1) x gdim
- tetrahedron_0 – (
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection_tetrahedron_tetrahedron
(const std::vector<Point> &tet_0, const std::vector<Point> &tet_1)¶ Parameters: - tet_0 –
- tet_1 –
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection_tetrahedron_triangle
(const MeshEntity &tetrahedron, const MeshEntity &triangle)¶ Compute triangulation of intersection of a tetrahedron and a triangle
Parameters: - tetrahedron – (
MeshEntity
) The tetrahedron. - triangle – (
MeshEntity
) The triangle
Returns: std::vector<double> A flattened array of simplices of dimension num_simplices x num_vertices x gdim = num_simplices x (tdim + 1) x gdim
- tetrahedron – (
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection_tetrahedron_triangle
(const std::vector<Point> &tet, const std::vector<Point> &tri)¶ Parameters: - tet –
- tri –
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection_triangle_interval
(const MeshEntity &triangle, const MeshEntity &interval)¶ Compute triangulation of intersection of a triangle and an interval
Parameters: - triangle – (
MeshEntity
) The triangle. - interval – (
MeshEntity
) The interval.
Returns: std::vector<double> A flattened array of simplices of dimension num_simplices x num_vertices x gdim = num_simplices x (tdim + 1) x gdim
- triangle – (
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection_triangle_interval
(const std::vector<Point> &triangle, const std::vector<Point> &interval, std::size_t gdim)¶ Parameters: - triangle –
- interval –
- gdim –
-
std::vector<double>
dolfin::IntersectionTriangulation::
triangulate_intersection_triangle_triangle
(const MeshEntity &triangle_0, const MeshEntity &triangle_1)¶ Compute triangulation of intersection of two triangles
Parameters: - triangle_0 – (
MeshEntity
) The first triangle. - triangle_1 – (
MeshEntity
) The second triangle.
Returns: std::vector<double> A flattened array of simplices of dimension num_simplices x num_vertices x gdim = num_simplices x (tdim + 1) x gdim
- triangle_0 – (
-
bool
MeshPointIntersection¶
C++ documentation for MeshPointIntersection
from dolfin/geometry/MeshPointIntersection.h
:
-
class
dolfin::
MeshPointIntersection
¶ This class represents an intersection between a
Mesh
and aPoint
. The resulting intersection is stored as a list of zero or more cells.-
dolfin::MeshPointIntersection::
MeshPointIntersection
(const Mesh &mesh, const Point &point)¶ Compute intersection between mesh and point.
Parameters: - mesh –
- point –
-
const std::vector<unsigned int> &
dolfin::MeshPointIntersection::
intersected_cells
() const¶ Return the list of (local) indices for intersected cells.
-
dolfin::MeshPointIntersection::
~MeshPointIntersection
()¶ Destructor.
-
Point¶
C++ documentation for Point
from dolfin/geometry/Point.h
:
-
class
dolfin::
Point
¶ A
Point
represents a point in \(\mathbb{R}^3\) with coordinates \(x, y, z,\) or alternatively, a vector in \(\mathbb{R}^3\) , supporting standard operations like the norm, distances, scalar and vector products etc.-
dolfin::Point::
Point
(const Array<double> &x)¶ Create point from
Array
Parameters: x – (Array<double>) Array
of coordinates.
-
dolfin::Point::
Point
(const Point &p)¶ Copy constructor
Parameters: p – ( Point()
) The object to be copied.
-
dolfin::Point::
Point
(const double x = 0.0, const double y = 0.0, const double z = 0.0)¶ Create a point at (x, y, z). Default value (0, 0, 0).
Parameters: - x – (double) The x-coordinate.
- y – (double) The y-coordinate.
- z – (double) The z-coordinate.
-
dolfin::Point::
Point
(std::size_t dim, const double *x)¶ Create point from array
Parameters: - dim – (std::size_t) Dimension of the array.
- x – (double) The array to create a
Point()
from.
-
std::array<double, 3>
dolfin::Point::
array
() const¶ Return copy of coordinate array Returns list of double The coordinates.
-
double *
dolfin::Point::
coordinates
()¶ Return coordinate array
Returns: double* The coordinates.
-
const double *
dolfin::Point::
coordinates
() const¶ Return coordinate array (const. version)
Returns: double* The coordinates.
-
const Point
dolfin::Point::
cross
(const Point &p) const¶ Compute cross product with given vector
Parameters: p – ( Point
) Another point.Returns: Point
The cross product.
-
double
dolfin::Point::
distance
(const Point &p) const¶ Compute distance to given point
Point p1(0, 4, 0); Point p2(2, 0, 4); info("%g", p1.distance(p2));
Parameters: p – ( Point
) The point to compute distance to.Returns: double The distance.
-
double
dolfin::Point::
dot
(const Point &p) const¶ Compute dot product with given vector
Point p1(1.0, 4.0, 8.0); Point p2(2.0, 0.0, 0.0); info("%g", p1.dot(p2));
Parameters: p – ( Point
) Another point.Returns: double The dot product.
-
double
dolfin::Point::
norm
() const¶ Compute norm of point representing a vector from the origin
Point p(1.0, 2.0, 2.0); info("%g", p.norm());
Returns: double The (Euclidean) norm of the vector from the origin to the point.
-
Point
dolfin::Point::
operator*
(double a) const¶ Multiplication with scalar.
Parameters: a –
-
const Point &
dolfin::Point::
operator*=
(double a)¶ Incremental multiplication with scalar.
Parameters: a –
-
Point
dolfin::Point::
operator+
(const Point &p) const¶ Compute sum of two points
Parameters: p – ( Point
)Returns: Point
-
const Point &
dolfin::Point::
operator+=
(const Point &p)¶ Add given point.
Parameters: p –
-
Point
dolfin::Point::
operator-
()¶ Unary minus.
-
Point
dolfin::Point::
operator-
(const Point &p) const¶ Compute difference of two points
Parameters: p – ( Point
)Returns: Point
-
const Point &
dolfin::Point::
operator-=
(const Point &p)¶ Subtract given point.
Parameters: p –
-
Point
dolfin::Point::
operator/
(double a) const¶ Division by scalar.
Parameters: a –
-
const Point &
dolfin::Point::
operator/=
(double a)¶ Incremental division by scalar.
Parameters: a –
-
const Point &
dolfin::Point::
operator=
(const Point &p)¶ Assignment operator.
Parameters: p –
-
double &
dolfin::Point::
operator[]
(std::size_t i)¶ Return address of coordinate in direction i Returns
Parameters: i – (std::size_t) Direction. Returns: double Address of coordinate in the given direction.
-
double
dolfin::Point::
operator[]
(std::size_t i) const¶ Return coordinate in direction i
Parameters: i – (std::size_t) Direction. Returns: double The coordinate in the given direction.
-
Point
dolfin::Point::
rotate
(const Point &a, double theta) const¶ Rotate around a given axis
Parameters: - a – (
Point
) The axis to rotate around. Must be unit length. - theta – (double) The rotation angle.
Returns: Point
The rotated point.- a – (
-
double
dolfin::Point::
squared_distance
(const Point &p) const¶ Compute squared distance to given point
Parameters: p – ( Point
) The point to compute distance to.Returns: double The squared distance.
-
double
dolfin::Point::
squared_norm
() const¶ Compute norm of point representing a vector from the origin
Point p(1.0, 2.0, 2.0); info("%g", p.squared_norm());
Returns: double The squared (Euclidean) norm of the vector from the origin of the point.
-
std::string
dolfin::Point::
str
(bool verbose = false) const¶ Return informal string representation (pretty-print)
Parameters: verbose – (bool) Flag to turn on additional output. Returns: std::string An informal representation of the function space.
-
double
dolfin::Point::
x
() const¶ Return x-coordinate
Returns: double The x-coordinate.
-
double
dolfin::Point::
y
() const¶ Return y-coordinate
Returns: double The y-coordinate.
-
double
dolfin::Point::
z
() const¶ Return z-coordinate
Returns: double The z-coordinate.
-
dolfin::Point::
~Point
()¶ Destructor.
-
SimplexQuadrature¶
C++ documentation for SimplexQuadrature
from dolfin/geometry/SimplexQuadrature.h
:
-
class
dolfin::
SimplexQuadrature
¶ Quadrature on simplices.
-
std::pair<std::vector<double>, std::vector<double>>
dolfin::SimplexQuadrature::
compute_quadrature_rule
(const Cell &cell, std::size_t order)¶ Compute quadrature rule for cell. Arguments cell (
Cell
) The cell. order (std::size_t) The order of convergence of the quadrature rule. Returns std::pair<std::vector<double>, std::vector<double> > A flattened array of quadrature points and a corresponding array of quadrature weights.Parameters: - cell –
- order –
-
std::pair<std::vector<double>, std::vector<double>>
dolfin::SimplexQuadrature::
compute_quadrature_rule
(const double *coordinates, std::size_t tdim, std::size_t gdim, std::size_t order)¶ Compute quadrature rule for simplex. Arguments coordinates (double *) A flattened array of simplex coordinates of dimension num_vertices x gdim = (tdim + 1)*gdim. tdim (std::size_t) The topological dimension of the simplex. gdim (std::size_t) The geometric dimension. order (std::size_t) The order of convergence of the quadrature rule. Returns std::pair<std::vector<double>, std::vector<double> > A flattened array of quadrature points and a corresponding array of quadrature weights.
Parameters: - coordinates –
- tdim –
- gdim –
- order –
-
std::pair<std::vector<double>, std::vector<double>>
dolfin::SimplexQuadrature::
compute_quadrature_rule_interval
(const double *coordinates, std::size_t gdim, std::size_t order)¶ Compute quadrature rule for interval. Arguments coordinates (double *) A flattened array of simplex coordinates of dimension num_vertices x gdim = 2*gdim. gdim (std::size_t) The geometric dimension. order (std::size_t) The order of convergence of the quadrature rule. Returns std::pair<std::vector<double>, std::vector<double> > A flattened array of quadrature points and a corresponding array of quadrature weights.
Parameters: - coordinates –
- gdim –
- order –
-
std::pair<std::vector<double>, std::vector<double>>
dolfin::SimplexQuadrature::
compute_quadrature_rule_tetrahedron
(const double *coordinates, std::size_t gdim, std::size_t order)¶ Compute quadrature rule for tetrahedron. Arguments coordinates (double *) A flattened array of simplex coordinates of dimension num_vertices x gdim = 4*gdim. gdim (std::size_t) The geometric dimension. order (std::size_t) The order of convergence of the quadrature rule. Returns std::pair<std::vector<double>, std::vector<double> > A flattened array of quadrature points and a corresponding array of quadrature weights.
Parameters: - coordinates –
- gdim –
- order –
-
std::pair<std::vector<double>, std::vector<double>>
dolfin::SimplexQuadrature::
compute_quadrature_rule_triangle
(const double *coordinates, std::size_t gdim, std::size_t order)¶ Compute quadrature rule for triangle. Arguments coordinates (double *) A flattened array of simplex coordinates of dimension num_vertices x gdim = 3*gdim. gdim (std::size_t) The geometric dimension. order (std::size_t) The order of convergence of the quadrature rule. Returns std::pair<std::vector<double>, std::vector<double> > A flattened array of quadrature points and a corresponding array of quadrature weights.
Parameters: - coordinates –
- gdim –
- order –
-
std::pair<std::vector<double>, std::vector<double>>