CGAL 5.4.3 - 2D Generalized Barycentric Coordinates

Free functions to compute barycentric weights and coordinates.

Functions

template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::discrete_harmonic_weights_2 (const PointRange &polygon, const typename GeomTraits::Point_2 &query, OutIterator w_begin, const GeomTraits &traits, const Computation_policy_2 policy=Computation_policy_2::FAST_WITH_EDGE_CASES)
 computes 2D discrete harmonic weights. More...
 
template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::discrete_harmonic_coordinates_2 (const PointRange &polygon, const typename GeomTraits::Point_2 &query, OutIterator c_begin, const GeomTraits &traits, const Computation_policy_2 policy=Computation_policy_2::PRECISE_WITH_EDGE_CASES)
 computes 2D discrete harmonic coordinates. More...
 
template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::mean_value_weights_2 (const PointRange &polygon, const typename GeomTraits::Point_2 &query, OutIterator w_begin, const GeomTraits &traits, const Computation_policy_2 policy=Computation_policy_2::FAST_WITH_EDGE_CASES)
 computes 2D mean value weights. More...
 
template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::mean_value_coordinates_2 (const PointRange &polygon, const typename GeomTraits::Point_2 &query, OutIterator c_begin, const GeomTraits &traits, const Computation_policy_2 policy=Computation_policy_2::PRECISE_WITH_EDGE_CASES)
 computes 2D mean value coordinates. More...
 
template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::wachspress_weights_2 (const PointRange &polygon, const typename GeomTraits::Point_2 &query, OutIterator w_begin, const GeomTraits &traits, const Computation_policy_2 policy=Computation_policy_2::FAST_WITH_EDGE_CASES)
 computes 2D Wachspress weights. More...
 
template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::wachspress_coordinates_2 (const PointRange &polygon, const typename GeomTraits::Point_2 &query, OutIterator c_begin, const GeomTraits &traits, const Computation_policy_2 policy=Computation_policy_2::PRECISE_WITH_EDGE_CASES)
 computes 2D Wachspress coordinates. More...
 
template<typename VertexRange , typename OutIterator , typename GeomTraits , typename PointMap >
std::pair< OutIterator, bool > CGAL::Barycentric_coordinates::boundary_coordinates_2 (const VertexRange &polygon, const typename GeomTraits::Point_2 &query, OutIterator c_begin, const GeomTraits &traits, const PointMap point_map)
 computes 2D boundary coordinates. More...
 
template<typename VertexRange , typename Query_2 , typename OutIterator , typename PointMap = CGAL::Identity_property_map<Query_2>>
std::pair< OutIterator, bool > CGAL::Barycentric_coordinates::boundary_coordinates_2 (const VertexRange &polygon, const Query_2 &query, OutIterator c_begin, const PointMap point_map=PointMap())
 computes 2D boundary coordinates. More...
 
template<typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::segment_coordinates_2 (const typename GeomTraits::Point_2 &p0, const typename GeomTraits::Point_2 &p1, const typename GeomTraits::Point_2 &query, OutIterator c_begin, const GeomTraits &traits)
 computes segment coordinates. More...
 
template<typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::triangle_coordinates_2 (const typename GeomTraits::Point_2 &p0, const typename GeomTraits::Point_2 &p1, const typename GeomTraits::Point_2 &p2, const typename GeomTraits::Point_2 &query, OutIterator c_begin, const GeomTraits &traits)
 computes triangle coordinates. More...
 

Function Documentation

◆ boundary_coordinates_2() [1/2]

template<typename VertexRange , typename OutIterator , typename GeomTraits , typename PointMap >
std::pair<OutIterator, bool> CGAL::Barycentric_coordinates::boundary_coordinates_2 ( const VertexRange &  polygon,
const typename GeomTraits::Point_2 &  query,
OutIterator  c_begin,
const GeomTraits &  traits,
const PointMap  point_map 
)

#include <CGAL/Barycentric_coordinates_2/boundary_coordinates_2.h>

computes 2D boundary coordinates.

This function computes boundary barycentric coordinates at a given query point with respect to the vertices of a simple polygon, that is one coordinate per vertex. The coordinates are stored in a destination range beginning at c_begin.

If query is at the vertex, the corresponding coordinate is set to one, while all other coordinates are zero. If query is on the edge, the two corresponding coordinates are segment coordinates, while all other coordinates are set to zero. If query is not on the boundary, all the coordinates are set to zero.

Internally, segment_coordinates_2() are used.

Template Parameters
VertexRangea model of ConstRange whose iterator type is RandomAccessIterator
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
GeomTraitsa model of BarycentricTraits_2
PointMapa model of ReadablePropertyMap whose key type is VertexRange::value_type and value type is GeomTraits::Point_2
Parameters
polygonan instance of VertexRange with 2D points, which form a simple polygon
querya query point
c_beginthe beginning of the destination range with the computed coordinates
traitsa traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type
point_mapan instance of PointMap that maps a vertex from polygon to Point_2
Returns
an output iterator to the element in the destination range, one past the last coordinate stored + the flag indicating whether the query point belongs to the polygon boundary
Precondition
polygon.size() >= 3
Examples:
Barycentric_coordinates_2/discrete_harmonic_coordinates.cpp.

◆ boundary_coordinates_2() [2/2]

template<typename VertexRange , typename Query_2 , typename OutIterator , typename PointMap = CGAL::Identity_property_map<Query_2>>
std::pair<OutIterator, bool> CGAL::Barycentric_coordinates::boundary_coordinates_2 ( const VertexRange &  polygon,
const Query_2 &  query,
OutIterator  c_begin,
const PointMap  point_map = PointMap() 
)

#include <CGAL/Barycentric_coordinates_2/boundary_coordinates_2.h>

computes 2D boundary coordinates.

This is an overload of the function boundary_coordinates_2.

Template Parameters
VertexRangea model of ConstRange whose iterator type is RandomAccessIterator
Query_2a model of Kernel::Point_2
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
PointMapa model of ReadablePropertyMap whose key type is VertexRange::value_type and value type is Query_2. The default is CGAL::Identity_property_map.
Parameters
polygonan instance of VertexRange with 2D points, which form a simple polygon
querya query point
c_beginthe beginning of the destination range with the computed coordinates
point_mapan instance of PointMap that maps a vertex from polygon to Query_2; the default initialization is provided
Returns
an output iterator to the element in the destination range, one past the last coordinate stored + the flag indicating whether the query point belongs to the polygon boundary
Precondition
polygon.size() >= 3

◆ discrete_harmonic_coordinates_2()

template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::discrete_harmonic_coordinates_2 ( const PointRange &  polygon,
const typename GeomTraits::Point_2 &  query,
OutIterator  c_begin,
const GeomTraits &  traits,
const Computation_policy_2  policy = Computation_policy_2::PRECISE_WITH_EDGE_CASES 
)

#include <CGAL/Barycentric_coordinates_2/Discrete_harmonic_coordinates_2.h>

computes 2D discrete harmonic coordinates.

This function computes 2D discrete harmonic coordinates at a given query point with respect to the vertices of a strictly convex polygon, that is one coordinate per vertex. The coordinates are stored in a destination range beginning at c_begin.

Internally, the class Discrete_harmonic_coordinates_2 is used. If one wants to process multiple query points, it is better to use that class. When using the free function, internal memory is allocated for each query point, while when using the class, it is allocated only once, which is much more efficient. However, for a few query points, it is easier to use this function. It can also be used when the processing time is not a concern.

Template Parameters
PointRangea model of ConstRange whose iterator type is RandomAccessIterator and value type is GeomTraits::Point_2
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
GeomTraitsa model of BarycentricTraits_2
Parameters
polygonan instance of PointRange with 2D points, which form a strictly convex polygon
querya query point
c_beginthe beginning of the destination range with the computed coordinates
traitsa traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type
policyone of the Computation_policy_2; the default is Computation_policy_2::PRECISE_WITH_EDGE_CASES
Returns
an output iterator to the element in the destination range, one past the last coordinate stored
Precondition
polygon.size() >= 3
polygon is simple
polygon is strictly convex

◆ discrete_harmonic_weights_2()

template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::discrete_harmonic_weights_2 ( const PointRange &  polygon,
const typename GeomTraits::Point_2 &  query,
OutIterator  w_begin,
const GeomTraits &  traits,
const Computation_policy_2  policy = Computation_policy_2::FAST_WITH_EDGE_CASES 
)

#include <CGAL/Barycentric_coordinates_2/Discrete_harmonic_coordinates_2.h>

computes 2D discrete harmonic weights.

This function computes 2D discrete harmonic weights at a given query point with respect to the vertices of a strictly convex polygon, that is one weight per vertex. The weights are stored in a destination range beginning at w_begin.

Internally, the class Discrete_harmonic_coordinates_2 is used. If one wants to process multiple query points, it is better to use that class. When using the free function, internal memory is allocated for each query point, while when using the class, it is allocated only once, which is much more efficient. However, for a few query points, it is easier to use this function. It can also be used when the processing time is not a concern.

Template Parameters
PointRangea model of ConstRange whose iterator type is RandomAccessIterator and value type is GeomTraits::Point_2
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
GeomTraitsa model of BarycentricTraits_2
Parameters
polygonan instance of PointRange with 2D points, which form a strictly convex polygon
querya query point
w_beginthe beginning of the destination range with the computed weights
traitsa traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type
policyone of the Computation_policy_2; the default is Computation_policy_2::FAST_WITH_EDGE_CASES
Returns
an output iterator to the element in the destination range, one past the last weight stored
Precondition
polygon.size() >= 3
polygon is simple
polygon is strictly convex

◆ mean_value_coordinates_2()

template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::mean_value_coordinates_2 ( const PointRange &  polygon,
const typename GeomTraits::Point_2 &  query,
OutIterator  c_begin,
const GeomTraits &  traits,
const Computation_policy_2  policy = Computation_policy_2::PRECISE_WITH_EDGE_CASES 
)

#include <CGAL/Barycentric_coordinates_2/Mean_value_coordinates_2.h>

computes 2D mean value coordinates.

This function computes 2D mean value coordinates at a given query point with respect to the vertices of a simple polygon, that is one coordinate per vertex. The coordinates are stored in a destination range beginning at c_begin.

Internally, the class Mean_value_coordinates_2 is used. If one wants to process multiple query points, it is better to use that class. When using the free function, internal memory is allocated for each query point, while when using the class, it is allocated only once, which is much more efficient. However, for a few query points, it is easier to use this function. It can also be used when the processing time is not a concern.

Template Parameters
PointRangea model of ConstRange whose iterator type is RandomAccessIterator and value type is GeomTraits::Point_2
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
GeomTraitsa model of BarycentricTraits_2
Parameters
polygonan instance of PointRange with 2D points, which form a simple polygon
querya query point
c_beginthe beginning of the destination range with the computed coordinates
traitsa traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type
policyone of the Computation_policy_2; the default is Computation_policy_2::PRECISE_WITH_EDGE_CASES
Returns
an output iterator to the element in the destination range, one past the last coordinate stored
Precondition
polygon.size() >= 3
polygon is simple
Examples:
Barycentric_coordinates_2/mean_value_coordinates.cpp, and Barycentric_coordinates_2/terrain_height_modeling.cpp.

◆ mean_value_weights_2()

template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::mean_value_weights_2 ( const PointRange &  polygon,
const typename GeomTraits::Point_2 &  query,
OutIterator  w_begin,
const GeomTraits &  traits,
const Computation_policy_2  policy = Computation_policy_2::FAST_WITH_EDGE_CASES 
)

#include <CGAL/Barycentric_coordinates_2/Mean_value_coordinates_2.h>

computes 2D mean value weights.

This function computes 2D mean value weights at a given query point with respect to the vertices of a simple polygon, that is one weight per vertex. The weights are stored in a destination range beginning at w_begin.

Internally, the class Mean_value_coordinates_2 is used. If one wants to process multiple query points, it is better to use that class. When using the free function, internal memory is allocated for each query point, while when using the class, it is allocated only once, which is much more efficient. However, for a few query points, it is easier to use this function. It can also be used when the processing time is not a concern.

Template Parameters
PointRangea model of ConstRange whose iterator type is RandomAccessIterator and value type is GeomTraits::Point_2
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
GeomTraitsa model of BarycentricTraits_2
Parameters
polygonan instance of PointRange with 2D points, which form a simple polygon
querya query point
w_beginthe beginning of the destination range with the computed weights
traitsa traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type
policyone of the Computation_policy_2; the default is Computation_policy_2::FAST_WITH_EDGE_CASES
Returns
an output iterator to the element in the destination range, one past the last weight stored
Precondition
polygon.size() >= 3
polygon is simple

◆ segment_coordinates_2()

template<typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::segment_coordinates_2 ( const typename GeomTraits::Point_2 &  p0,
const typename GeomTraits::Point_2 &  p1,
const typename GeomTraits::Point_2 &  query,
OutIterator  c_begin,
const GeomTraits &  traits 
)

#include <CGAL/Barycentric_coordinates_2/segment_coordinates_2.h>

computes segment coordinates.

This function computes barycentric coordinates at a given query point with respect to the end points p0 and p1 of a segment that is one coordinate per end point. The coordinates are stored in a destination range beginning at c_begin.

After the coordinates \(b_0\) and \(b_1\) are computed, the query point \(q\) can be obtained as \(q = b_0p_0 + b_1p_1\). If \(q\) does not belong to the line through \(p_0\) and \(p_1\), it is projected onto this line, and only then the coordinates are computed. See more details in the user manual here.

Template Parameters
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
GeomTraitsa model of BarycentricTraits_2
Parameters
p0the first vertex of a segment
p1the second vertex of a segment
querya query point
c_beginthe beginning of the destination range with the computed coordinates
traitsa traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type
Returns
an output iterator to the element in the destination range, one past the last coordinate stored
Precondition
p0 != p1
Examples:
Barycentric_coordinates_2/segment_coordinates.cpp.

◆ triangle_coordinates_2()

template<typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::triangle_coordinates_2 ( const typename GeomTraits::Point_2 &  p0,
const typename GeomTraits::Point_2 &  p1,
const typename GeomTraits::Point_2 &  p2,
const typename GeomTraits::Point_2 &  query,
OutIterator  c_begin,
const GeomTraits &  traits 
)

#include <CGAL/Barycentric_coordinates_2/triangle_coordinates_2.h>

computes triangle coordinates.

This function computes barycentric coordinates at a given query point with respect to the points p0, p1, and p2, which form a triangle, that is one coordinate per point. The coordinates are stored in a destination range beginning at c_begin.

After the coordinates \(b_0\), \(b_1\), and \(b_2\) are computed, the query point \(q\) can be obtained as \(q = b_0p_0 + b_1p_1 + b_2p_2\). See more details in the user manual here.

Template Parameters
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
GeomTraitsa model of BarycentricTraits_2
Parameters
p0the first vertex of a triangle
p1the second vertex of a triangle
p2the third vertex of a triangle
querya query point
c_beginthe beginning of the destination range with the computed coordinates
traitsa traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type
Returns
an output iterator to the element in the destination range, one past the last coordinate stored
Precondition
area_2(p0, p1, p2) != 0
Examples:
Barycentric_coordinates_2/triangle_coordinates.cpp.

◆ wachspress_coordinates_2()

template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::wachspress_coordinates_2 ( const PointRange &  polygon,
const typename GeomTraits::Point_2 &  query,
OutIterator  c_begin,
const GeomTraits &  traits,
const Computation_policy_2  policy = Computation_policy_2::PRECISE_WITH_EDGE_CASES 
)

#include <CGAL/Barycentric_coordinates_2/Wachspress_coordinates_2.h>

computes 2D Wachspress coordinates.

This function computes 2D Wachspress coordinates at a given query point with respect to the vertices of a strictly convex polygon, that is one coordinate per vertex. The coordinates are stored in a destination range beginning at c_begin.

Internally, the class Wachspress_coordinates_2 is used. If one wants to process multiple query points, it is better to use that class. When using the free function, internal memory is allocated for each query point, while when using the class, it is allocated only once, which is much more efficient. However, for a few query points, it is easier to use this function. It can also be used when the processing time is not a concern.

Template Parameters
PointRangea model of ConstRange whose iterator type is RandomAccessIterator and value type is GeomTraits::Point_2
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
GeomTraitsa model of BarycentricTraits_2
Parameters
polygonan instance of PointRange with 2D points, which form a strictly convex polygon
querya query point
c_beginthe beginning of the destination range with the computed coordinates
traitsa traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type
policyone of the Computation_policy_2; the default is Computation_policy_2::PRECISE_WITH_EDGE_CASES
Returns
an output iterator to the element in the destination range, one past the last coordinate stored
Precondition
polygon.size() >= 3
polygon is simple
polygon is strictly convex
Examples:
Barycentric_coordinates_2/wachspress_coordinates.cpp.

◆ wachspress_weights_2()

template<typename PointRange , typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::wachspress_weights_2 ( const PointRange &  polygon,
const typename GeomTraits::Point_2 &  query,
OutIterator  w_begin,
const GeomTraits &  traits,
const Computation_policy_2  policy = Computation_policy_2::FAST_WITH_EDGE_CASES 
)

#include <CGAL/Barycentric_coordinates_2/Wachspress_coordinates_2.h>

computes 2D Wachspress weights.

This function computes 2D Wachspress weights at a given query point with respect to the vertices of a strictly convex polygon, that is one weight per vertex. The weights are stored in a destination range beginning at w_begin.

Internally, the class Wachspress_coordinates_2 is used. If one wants to process multiple query points, it is better to use that class. When using the free function, internal memory is allocated for each query point, while when using the class, it is allocated only once, which is much more efficient. However, for a few query points, it is easier to use this function. It can also be used when the processing time is not a concern.

Template Parameters
PointRangea model of ConstRange whose iterator type is RandomAccessIterator and value type is GeomTraits::Point_2
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
GeomTraitsa model of BarycentricTraits_2
Parameters
polygonan instance of PointRange with 2D points, which form a strictly convex polygon
querya query point
w_beginthe beginning of the destination range with the computed weights
traitsa traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type
policyone of the Computation_policy_2; the default is Computation_policy_2::FAST_WITH_EDGE_CASES
Returns
an output iterator to the element in the destination range, one past the last weight stored
Precondition
polygon.size() >= 3
polygon is simple
polygon is strictly convex