dolfin/geometry

Documentation for C++ code found in dolfin/geometry/*.h

Functions

exactinit

C++ documentation for exactinit from dolfin/geometry/predicates.h:

void dolfin::exactinit()

Initialize tolerances for exact arithmetic.

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 and Point . 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:

Point dolfin::operator*(double a, const Point &p)

Multiplication with scalar.

Parameters:
  • a
  • p

operator<<

C++ documentation for operator<< from dolfin/geometry/Point.h:

std::ostream &dolfin::operator<<(std::ostream &stream, const Point &point)

Output of Point to stream.

Parameters:
  • stream
  • point

orient1d

C++ documentation for orient1d from dolfin/geometry/predicates.h:

double dolfin::orient1d(double a, double b, double x)

Compute relative orientation of point x wrt segment [a, b].

Parameters:
  • a
  • b
  • x

orient2d

C++ documentation for orient2d from dolfin/geometry/predicates.h:

double dolfin::orient2d(const Point &a, const Point &b, const Point &c)

Convenience function using dolfin::Point .

Parameters:
  • a
  • b
  • c

orient3d

C++ documentation for orient3d from dolfin/geometry/predicates.h:

double dolfin::orient3d(const Point &a, const Point &b, const Point &c, const Point &d)

Convenience function using dolfin::Point .

Parameters:
  • a
  • b
  • c
  • d

Variables

predicate_initialization

C++ documentation for predicate_initialization from dolfin/geometry/predicates.cpp:

PredicateInitialization dolfin::predicate_initialization

Initialize the predicate.

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 entity i in the first list collides with entity i 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 entity i in the first list collides with entity i 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 the Point . Returns std::vector<unsigned int> A list of process numbers where the Mesh 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

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

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

CollisionPredicates

C++ documentation for CollisionPredicates from dolfin/geometry/CollisionPredicates.h:

class dolfin::CollisionPredicates

This class implements algorithms for detecting pairwise collisions between mesh entities of varying dimensions.

bool dolfin::CollisionPredicates::collides(const MeshEntity &entity, const Point &point)

Check whether entity collides with point. Arguments entity (MeshEntity ) The entity. point (Point ) The point. Returns bool True iff entity collides with cell.

Parameters:
  • entity
  • point
bool dolfin::CollisionPredicates::collides(const MeshEntity &entity_0, const MeshEntity &entity_1)

Check whether two entities collide. Arguments entity_0 (MeshEntity ) The first entity. entity_1 (MeshEntity ) The second entity. Returns bool True iff entity collides with cell.

Parameters:
  • entity_0
  • entity_1
bool dolfin::CollisionPredicates::collides_segment_point(const Point &p0, const Point &p1, const Point &point, std::size_t gdim)

Check whether segment p0-p1 collides with point.

Parameters:
  • p0
  • p1
  • point
  • gdim
bool dolfin::CollisionPredicates::collides_segment_point_1d(double p0, double p1, double point)

Check whether segment p0-p1 collides with point (1D version)

Parameters:
  • p0
  • p1
  • point
bool dolfin::CollisionPredicates::collides_segment_point_2d(const Point &p0, const Point &p1, const Point &point)

Check whether segment p0-p1 collides with point (2D version)

Parameters:
  • p0
  • p1
  • point
bool dolfin::CollisionPredicates::collides_segment_point_3d(const Point &p0, const Point &p1, const Point &point)

Check whether segment p0-p1 collides with point (3D version)

Parameters:
  • p0
  • p1
  • point
bool dolfin::CollisionPredicates::collides_segment_segment(const Point &p0, const Point &p1, const Point &q0, const Point &q1, std::size_t gdim)

Check whether segment p0-p1 collides with segment q0-q1.

Parameters:
  • p0
  • p1
  • q0
  • q1
  • gdim
bool dolfin::CollisionPredicates::collides_segment_segment_1d(double p0, double p1, double q0, double q1)

Check whether segment p0-p1 collides with segment q0-q1 (1D version)

Parameters:
  • p0
  • p1
  • q0
  • q1
bool dolfin::CollisionPredicates::collides_segment_segment_2d(const Point &p0, const Point &p1, const Point &q0, const Point &q1)

Check whether segment p0-p1 collides with segment q0-q1 (2D version)

Parameters:
  • p0
  • p1
  • q0
  • q1
bool dolfin::CollisionPredicates::collides_segment_segment_3d(const Point &p0, const Point &p1, const Point &q0, const Point &q1)

Check whether segment p0-p1 collides with segment q0-q1 (3D version)

Parameters:
  • p0
  • p1
  • q0
  • q1
bool dolfin::CollisionPredicates::collides_tetrahedron_point_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &point)

Check whether tetrahedron p0-p1-p2-p3 collides with point.

Parameters:
  • p0
  • p1
  • p2
  • p3
  • point
bool dolfin::CollisionPredicates::collides_tetrahedron_segment_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1)

Check whether tetrahedron p0-p1-p2-p3 collides with segment q0-q1.

Parameters:
  • p0
  • p1
  • p2
  • p3
  • q0
  • q1
bool dolfin::CollisionPredicates::collides_tetrahedron_tetrahedron_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2, const Point &q3)

Check whether tetrahedron p0-p1-p2-p3 collides with tetrahedron q0-q1-q2.

Parameters:
  • p0
  • p1
  • p2
  • p3
  • q0
  • q1
  • q2
  • q3
bool dolfin::CollisionPredicates::collides_tetrahedron_triangle_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2)

Check whether tetrahedron p0-p1-p2-p3 collides with triangle q0-q1-q2.

Parameters:
  • p0
  • p1
  • p2
  • p3
  • q0
  • q1
  • q2
bool dolfin::CollisionPredicates::collides_triangle_point(const Point &p0, const Point &p1, const Point &p2, const Point &point, std::size_t gdim)

Check whether triangle p0-p1-p2 collides with point.

Parameters:
  • p0
  • p1
  • p2
  • point
  • gdim
bool dolfin::CollisionPredicates::collides_triangle_point_2d(const Point &p0, const Point &p1, const Point &p2, const Point &point)

Check whether triangle p0-p1-p2 collides with point (2D version)

Parameters:
  • p0
  • p1
  • p2
  • point
bool dolfin::CollisionPredicates::collides_triangle_point_3d(const Point &p0, const Point &p1, const Point &p2, const Point &point)

Check whether triangle p0-p1-p2 collides with point (3D version)

Parameters:
  • p0
  • p1
  • p2
  • point
bool dolfin::CollisionPredicates::collides_triangle_segment(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, std::size_t gdim)

Check whether triangle p0-p1-p2 collides with segment q0-q1.

Parameters:
  • p0
  • p1
  • p2
  • q0
  • q1
  • gdim
bool dolfin::CollisionPredicates::collides_triangle_segment_2d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)

Check whether triangle p0-p1-p2 collides with segment q0-q1 (2D version)

Parameters:
  • p0
  • p1
  • p2
  • q0
  • q1
bool dolfin::CollisionPredicates::collides_triangle_segment_3d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)

Check whether triangle p0-p1-p2 collides with segment q0-q1 (3D version)

Parameters:
  • p0
  • p1
  • p2
  • q0
  • q1
bool dolfin::CollisionPredicates::collides_triangle_triangle(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2, std::size_t gdim)

Check whether triangle p0-p1-p2 collides with triangle q0-q1-q2.

Parameters:
  • p0
  • p1
  • p2
  • q0
  • q1
  • q2
  • gdim
bool dolfin::CollisionPredicates::collides_triangle_triangle_2d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)

Check whether triangle p0-p1-p2 collides with triangle q0-q1-q2 (2D version)

Parameters:
  • p0
  • p1
  • p2
  • q0
  • q1
  • q2
bool dolfin::CollisionPredicates::collides_triangle_triangle_3d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)

Check whether triangle p0-p1-p2 collides with triangle q0-q1-q2 (3D version)

Parameters:
  • p0
  • p1
  • p2
  • q0
  • q1
  • q2

ConvexTriangulation

C++ documentation for ConvexTriangulation from dolfin/geometry/ConvexTriangulation.h:

class dolfin::ConvexTriangulation

This class implements algorithms for triangulating convex domains represented as a set of points.

bool dolfin::ConvexTriangulation::selfintersects(const std::vector<std::vector<Point>> &p)

Determine if there are self-intersecting tetrahedra.

Parameters:p
std::vector<std::vector<Point>> dolfin::ConvexTriangulation::triangulate(const std::vector<Point> &p, std::size_t gdim, std::size_t tdim)

Tdim independent wrapper.

Parameters:
  • p
  • gdim
  • tdim
std::vector<std::vector<Point>> dolfin::ConvexTriangulation::triangulate_1d(const std::vector<Point> &pm, std::size_t gdim)

Triangulate 1D.

Parameters:
  • pm
  • gdim
std::vector<std::vector<Point>> dolfin::ConvexTriangulation::triangulate_graham_scan_2d(const std::vector<Point> &pm)

Triangulate using the Graham scan 2D.

Parameters:pm
std::vector<std::vector<Point>> dolfin::ConvexTriangulation::triangulate_graham_scan_3d(const std::vector<Point> &pm)

Triangulate using the Graham scan 3D.

Parameters:pm

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.

GeometryDebugging

C++ documentation for GeometryDebugging from dolfin/geometry/GeometryDebugging.h:

class dolfin::GeometryDebugging

This class provides useful functionality for debugging algorithms dealing with geometry such as collision detection and intersection triangulation.

void dolfin::GeometryDebugging::init_plot()

Initialize plotting (print matplotlib code).

void dolfin::GeometryDebugging::plot(const Point &point)

Plot a point (print matplotlib code). Example usage: plot(p0)

Parameters:point
void dolfin::GeometryDebugging::plot(const std::vector<Point> &simplex)

Plot a simplex (print matplotlib code). Example usage: plot({p0, p1, p2})

Parameters:simplex
void dolfin::GeometryDebugging::plot(const std::vector<Point> &simplex_0, const std::vector<Point> &simplex_1)

Plot a pair of simplices (print matplotlib code). Example usage: plot({p0, p1, p2}, {q0, q1})

Parameters:
  • simplex_0
  • simplex_1
std::string dolfin::GeometryDebugging::point2string(const Point &p)

Compact point to string conversion.

Parameters:p
void dolfin::GeometryDebugging::print(const Point &point)

Print coordinates of a point. Example usage: print(p0)

Parameters:point
void dolfin::GeometryDebugging::print(const std::vector<Point> &simplex)

Print coordinates of a simplex. Example usage: print({p0, p1, p2})

Parameters:simplex
void dolfin::GeometryDebugging::print(const std::vector<Point> &simplex_0, const std::vector<Point> &simplex_1)

Print coordinates of a pair of simplices. Example usage: print({p0, p1, p2}, {q0, q1})

Parameters:
  • simplex_0
  • simplex_1
std::string dolfin::GeometryDebugging::simplex2string(const std::vector<Point> &simplex)

Compact simplex to string conversion.

Parameters:simplex

GeometryPredicates

C++ documentation for GeometryPredicates from dolfin/geometry/GeometryPredicates.h:

class dolfin::GeometryPredicates

This class implements geometric predicates, i.e. function that return either true or false.

bool dolfin::GeometryPredicates::convex_hull_is_degenerate(const std::vector<Point> &p, std::size_t gdim)

Check whether the convex hull is degenerate.

Parameters:
  • p
  • gdim
bool dolfin::GeometryPredicates::is_degenerate(const std::vector<Point> &simplex, std::size_t gdim)

Check whether simplex is degenerate.

Parameters:
  • simplex
  • gdim
bool dolfin::GeometryPredicates::is_degenerate_2d(const std::vector<Point> &simplex)

Check whether simplex is degenerate (2D version)

Parameters:simplex
bool dolfin::GeometryPredicates::is_degenerate_3d(const std::vector<Point> &simplex)

Check whether simplex is degenerate (3D version)

Parameters:simplex
bool dolfin::GeometryPredicates::is_finite(const std::vector<Point> &simplex)

Check whether simplex is finite (not Inf or NaN)

Parameters:simplex
bool dolfin::GeometryPredicates::is_finite(const std::vector<double> &simplex)

Check whether simplex is finite (not Inf or NaN)

Parameters:simplex

GeometryTools

C++ documentation for GeometryTools from dolfin/geometry/GeometryTools.h:

class dolfin::GeometryTools

This class provides useful tools (functions) for computational geometry.

bool dolfin::GeometryTools::contains(double a, double b, double x)

Check whether x in [a, b].

Parameters:
  • a
  • b
  • x
bool dolfin::GeometryTools::contains_strict(double a, double b, double x)

Check whether x in (a, b)

Parameters:
  • a
  • b
  • x
Point dolfin::GeometryTools::cross_product(const Point &a, const Point &b, const Point &c)

Compute numerically stable cross product (a - c) x (b - c)

Parameters:
  • a
  • b
  • c
double dolfin::GeometryTools::determinant(const Point &ab, const Point &dc, const Point &ec)

Compute determinant of 3 x 3 matrix defined by vectors, ab, dc, ec.

Parameters:
  • ab
  • dc
  • ec
std::size_t dolfin::GeometryTools::major_axis_2d(const Point &v)

Compute major (largest) axis of vector (2D)

Parameters:v
std::size_t dolfin::GeometryTools::major_axis_3d(const Point &v)

Compute major (largest) axis of vector (3D)

Parameters:v
double dolfin::GeometryTools::project_to_axis_2d(const Point &p, std::size_t axis)

Project point to axis (2D)

Parameters:
  • p
  • axis
Point dolfin::GeometryTools::project_to_plane_3d(const Point &p, std::size_t axis)

Project point to plane (3D)

Parameters:
  • p
  • axis

IntersectionConstruction

C++ documentation for IntersectionConstruction from dolfin/geometry/IntersectionConstruction.h:

class dolfin::IntersectionConstruction

This class implements algorithms for computing pairwise intersections of simplices. The computed intersection is always convex and represented as a set of points s.t. the intersection is the convex hull of these points.

std::vector<Point> dolfin::IntersectionConstruction::intersection(const MeshEntity &entity_0, const MeshEntity &entity_1)

Compute intersection of two entities. Arguments entity_0 (MeshEntity ) The first entity. entity_1 (MeshEntity ) The second entity. Returns std::vector<Pointdouble> A vector of points s.t. the intersection is the convex hull of these points.

Parameters:
  • entity_0
  • entity_1
std::vector<Point> dolfin::IntersectionConstruction::intersection(const std::vector<Point> &points_0, const std::vector<Point> &points_1, std::size_t gdim)

Compute intersection of two entities. This version takes two vectors of points representing the entities. Arguments points_0 (std::vector<Point>) The vertex coordinates of the first entity. points_1 (std::vector<Point>) The vertex coordinates of the second entity. gdim (std::size_t) The geometric dimension. Returns std::vector<Point> A vector of points s.t. the intersection is the convex hull of these points.

Parameters:
  • points_0
  • points_1
  • gdim
std::vector<double> dolfin::IntersectionConstruction::intersection_point_point_1d(double p0, double q0)

Compute intersection of points p0 and q0 (1D)

Parameters:
  • p0
  • q0
std::vector<Point> dolfin::IntersectionConstruction::intersection_point_point_2d(const Point &p0, const Point &q0)

Compute intersection of points p0 and q0 (2D)

Parameters:
  • p0
  • q0
std::vector<Point> dolfin::IntersectionConstruction::intersection_point_point_3d(const Point &p0, const Point &q0)

Compute intersection of points p0 and q0 (3D)

Parameters:
  • p0
  • q0
std::vector<double> dolfin::IntersectionConstruction::intersection_segment_point_1d(double p0, double p1, double q0)

Compute intersection of segment p0-p1 with point q0 (1D)

Parameters:
  • p0
  • p1
  • q0
std::vector<Point> dolfin::IntersectionConstruction::intersection_segment_point_2d(const Point &p0, const Point &p1, const Point &q0)

Compute intersection of segment p0-p1 with point q0 (2D)

Parameters:
  • p0
  • p1
  • q0
std::vector<Point> dolfin::IntersectionConstruction::intersection_segment_point_3d(const Point &p0, const Point &p1, const Point &q0)

Compute intersection of segment p0-p1 with point q0 (3D)

Parameters:
  • p0
  • p1
  • q0
std::vector<double> dolfin::IntersectionConstruction::intersection_segment_segment_1d(double p0, double p1, double q0, double q1)

Compute intersection of segment p0-p1 with segment q0-q1 (1D)

Parameters:
  • p0
  • p1
  • q0
  • q1
std::vector<Point> dolfin::IntersectionConstruction::intersection_segment_segment_2d(const Point &p0, const Point &p1, const Point &q0, const Point &q1)

Compute intersection of segment p0-p1 with segment q0-q1 (2D)

Parameters:
  • p0
  • p1
  • q0
  • q1
std::vector<Point> dolfin::IntersectionConstruction::intersection_segment_segment_3d(const Point &p0, const Point &p1, const Point &q0, const Point &q1)

Compute intersection of segment p0-p1 with segment q0-q1 (3D)

Parameters:
  • p0
  • p1
  • q0
  • q1
std::vector<Point> dolfin::IntersectionConstruction::intersection_tetrahedron_point_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0)

Compute intersection of tetrahedron p0-p1-p2-p3 with point q0 (3D)

Parameters:
  • p0
  • p1
  • p2
  • p3
  • q0
std::vector<Point> dolfin::IntersectionConstruction::intersection_tetrahedron_segment_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1)

Compute intersection of tetrahedron p0-p1-p2-p3 with segment q0-q1 (3D)

Parameters:
  • p0
  • p1
  • p2
  • p3
  • q0
  • q1
std::vector<Point> dolfin::IntersectionConstruction::intersection_tetrahedron_tetrahedron_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2, const Point &q3)

Compute intersection of tetrahedron p0-p1-p2-p3 with tetrahedron q0-q1-q2-q3 (3D)

Parameters:
  • p0
  • p1
  • p2
  • p3
  • q0
  • q1
  • q2
  • q3
std::vector<Point> dolfin::IntersectionConstruction::intersection_tetrahedron_triangle_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2)

Compute intersection of tetrahedron p0-p1-p2-p3 with triangle q0-q1-q2 (3D)

Parameters:
  • p0
  • p1
  • p2
  • p3
  • q0
  • q1
  • q2
std::vector<Point> dolfin::IntersectionConstruction::intersection_triangle_point_2d(const Point &p0, const Point &p1, const Point &p2, const Point &q0)

Compute intersection of triangle p0-p1-p2 with point q0 (2D)

Parameters:
  • p0
  • p1
  • p2
  • q0
std::vector<Point> dolfin::IntersectionConstruction::intersection_triangle_point_3d(const Point &p0, const Point &p1, const Point &p2, const Point &q0)

Compute intersection of triangle p0-p1-p2 with point q0 (3D)

Parameters:
  • p0
  • p1
  • p2
  • q0
std::vector<Point> dolfin::IntersectionConstruction::intersection_triangle_segment_2d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)

Compute intersection of triangle p0-p1-p2 with segment q0-q1 (2D)

Parameters:
  • p0
  • p1
  • p2
  • q0
  • q1
std::vector<Point> dolfin::IntersectionConstruction::intersection_triangle_segment_3d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)

Compute intersection of triangle p0-p1-p2 with segment q0-q1 (3D)

Parameters:
  • p0
  • p1
  • p2
  • q0
  • q1
std::vector<Point> dolfin::IntersectionConstruction::intersection_triangle_triangle_2d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)

Compute intersection of triangle p0-p1-p2 with triangle q0-q1-q2 (2D)

Parameters:
  • p0
  • p1
  • p2
  • q0
  • q1
  • q2
std::vector<Point> dolfin::IntersectionConstruction::intersection_triangle_triangle_3d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)

Compute intersection of triangle p0-p1-p2 with triangle q0-q1-q2 (3D)

Parameters:
  • p0
  • p1
  • p2
  • q0
  • q1
  • q2

MeshPointIntersection

C++ documentation for MeshPointIntersection from dolfin/geometry/MeshPointIntersection.h:

class dolfin::MeshPointIntersection

This class represents an intersection between a Mesh and a Point . 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.
bool dolfin::Point::operator!=(const Point &p) const

Not equal to operator.

Parameters:p
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
bool dolfin::Point::operator==(const Point &p) const

Equal to 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.

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.

PredicateInitialization

C++ documentation for PredicateInitialization from dolfin/geometry/predicates.h:

class dolfin::PredicateInitialization

Class used for automatic initialization of tolerances at startup. A global instance is defined inside predicates.cpp to ensure that the constructor and thus exactinit() is called.

dolfin::PredicateInitialization::PredicateInitialization()

SimplexQuadrature

C++ documentation for SimplexQuadrature from dolfin/geometry/SimplexQuadrature.h:

class dolfin::SimplexQuadrature

This class defines quadrature rules for simplices.

Eigen::MatrixXd dolfin::SimplexQuadrature::Chebyshev_Vandermonde_matrix(const std::pair<std::vector<double>, std::vector<double>> &qr, std::size_t gdim, std::size_t N)
Parameters:
  • qr
  • gdim
  • N
std::vector<Eigen::VectorXd> dolfin::SimplexQuadrature::Chebyshev_polynomial(const Eigen::VectorXd &x, std::size_t N)
Parameters:
  • x
  • N
dolfin::SimplexQuadrature::SimplexQuadrature(std::size_t tdim, std::size_t order)

Create SimplexQuadrature() rules for reference simplex Arguments tdim (std::size_t) The topological dimension of the simplex. order (std::size_t) The order of convergence of the quadrature rule.

Parameters:
  • tdim
  • order
std::size_t dolfin::SimplexQuadrature::choose(std::size_t n, std::size_t k)
Parameters:
  • n
  • k
std::vector<std::size_t> dolfin::SimplexQuadrature::compress(std::pair<std::vector<double>, std::vector<double>> &qr, std::size_t gdim, std::size_t quadrature_order)

Compress a quadrature rule using algorithms from Compression of multivariate discrete measures and applications A. Sommariva, M. Vianello Numerical Functional Analysis and Optimization Volume 36, 2015 - Issue 9 Arguments qr (std::pair<std::vector<double>, std::vector<double>>) The quadrature rule to be compressed gdim (std::size_t) The geometric dimension quadrature_order (std::size_t) The order of the quadrature rule Returns std::vector<std::size_t> The indices of the points that were kept (empty if no compression was made)

Parameters:
  • qr
  • gdim
  • quadrature_order
std::pair<std::vector<double>, std::vector<double>> dolfin::SimplexQuadrature::compute_quadrature_rule(const Cell &cell, std::size_t order) const

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 std::vector<Point> &coordinates, std::size_t gdim, std::size_t order) const

Compute quadrature rule for simplex. Arguments coordinates (std::vector<Point>) Vertex coordinates for the simplex 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
  • gdim
  • order
std::pair<std::vector<double>, std::vector<double>> dolfin::SimplexQuadrature::compute_quadrature_rule_interval(const std::vector<Point> &coordinates, std::size_t gdim, std::size_t order) const

Compute quadrature rule for interval. Arguments coordinates (std::vector<Point>) Vertex coordinates for 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
  • gdim
  • order
std::pair<std::vector<double>, std::vector<double>> dolfin::SimplexQuadrature::compute_quadrature_rule_tetrahedron(const std::vector<Point> &coordinates, std::size_t gdim, std::size_t order) const

Compute quadrature rule for tetrahedron. Arguments coordinates (std::vector<Point>) Vertex coordinates for 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
  • gdim
  • order
std::pair<std::vector<double>, std::vector<double>> dolfin::SimplexQuadrature::compute_quadrature_rule_triangle(const std::vector<Point> &coordinates, std::size_t gdim, std::size_t order) const

Compute quadrature rule for triangle. Arguments coordinates (std::vector<Point>) Vertex coordinates for 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
  • gdim
  • order
std::size_t dolfin::SimplexQuadrature::dunavant_order_num(std::size_t rule)
Parameters:rule
void dolfin::SimplexQuadrature::dunavant_rule(std::size_t order, std::vector<std::vector<double>> &p, std::vector<double> &w)
Parameters:
  • order
  • p
  • w
std::vector<std::size_t> dolfin::SimplexQuadrature::dunavant_suborder(int rule, int suborder_num)
Parameters:
  • rule
  • suborder_num
std::size_t dolfin::SimplexQuadrature::dunavant_suborder_num(int rule)
Parameters:rule
void dolfin::SimplexQuadrature::dunavant_subrule(std::size_t rule, std::size_t suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &w)
Parameters:
  • rule
  • suborder_num
  • suborder_xyz
  • w
void dolfin::SimplexQuadrature::dunavant_subrule_01(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_02(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_03(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_04(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_05(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_06(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_07(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_08(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_09(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_10(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_11(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_12(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_13(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_14(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_15(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_16(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_17(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_18(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_19(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
void dolfin::SimplexQuadrature::dunavant_subrule_20(int suborder_num, std::vector<double> &suborder_xyz, std::vector<double> &suborder_w)
Parameters:
  • suborder_num
  • suborder_xyz
  • suborder_w
std::vector<std::vector<std::size_t>> dolfin::SimplexQuadrature::grlex(std::size_t gdim, std::size_t N)
Parameters:
  • gdim
  • N
int dolfin::SimplexQuadrature::i4_modp(int i, int j)
Parameters:
  • i
  • j
int dolfin::SimplexQuadrature::i4_wrap(int ival, int ilo, int ihi)
Parameters:
  • ival
  • ilo
  • ihi
void dolfin::SimplexQuadrature::legendre_compute_glr(std::size_t n, std::vector<double> &x, std::vector<double> &w)
Parameters:
  • n
  • x
  • w
void dolfin::SimplexQuadrature::legendre_compute_glr0(std::size_t n, double &p, double &pp)
Parameters:
  • n
  • p
  • pp
void dolfin::SimplexQuadrature::legendre_compute_glr1(std::size_t n, std::vector<double> &x, std::vector<double> &w)
Parameters:
  • n
  • x
  • w
void dolfin::SimplexQuadrature::legendre_compute_glr2(double pn0, int n, double &x1, double &d1)
Parameters:
  • pn0
  • n
  • x1
  • d1
double dolfin::SimplexQuadrature::rk2_leg(double t1, double t2, double x, int n)
Parameters:
  • t1
  • t2
  • x
  • n
void dolfin::SimplexQuadrature::setup_qr_reference_interval(std::size_t order)
Parameters:order
void dolfin::SimplexQuadrature::setup_qr_reference_tetrahedron(std::size_t order)
Parameters:order
void dolfin::SimplexQuadrature::setup_qr_reference_triangle(std::size_t order)
Parameters:order
double dolfin::SimplexQuadrature::ts_mult(std::vector<double> &u, double h, int n)
Parameters:
  • u
  • h
  • n