CGAL 6.0 - Vector Graphics on Triangulated Surface Meshes
|
Functions | |
template<class FT , class TriangleMesh , class EdgeLocationRange > | |
void | CGAL::Vector_graphics_on_surfaces::locally_shortest_path (CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, FT > src, CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, FT > tgt, const TriangleMesh &tmesh, EdgeLocationRange &edge_locations, const Dual_geodesic_solver< FT > &solver) |
computes an approximated geodesic shortest path between two locations on a tmesh . | |
template<class TriangleMesh , class FT > | |
std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, FT > > | CGAL::Vector_graphics_on_surfaces::recursive_de_Casteljau (const TriangleMesh &mesh, const Bezier_segment< TriangleMesh, FT > &control_points, const int num_subdiv, const Dual_geodesic_solver< FT > &solver=Dual_geodesic_solver< FT >()) |
computes a discretization of a Bézier segment defined by the location of four control points on tmesh . | |
template<class K , class TriangleMesh > | |
std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > | CGAL::Vector_graphics_on_surfaces::straightest_geodesic (const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > &src, const typename K::Vector_2 &dir, const typename K::FT len, const TriangleMesh &tmesh) |
computes a path on a triangle mesh that is computed by starting a walk on tmesh given a direction and a maximum distance. | |
template<class K , class TriangleMesh > | |
std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > | CGAL::Vector_graphics_on_surfaces::trace_geodesic_polygon (const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > ¢er, const std::vector< typename K::Vector_2 > &directions, const std::vector< typename K::FT > &lengths, const TriangleMesh &tmesh, const Dual_geodesic_solver< typename K::FT > &solver={}) |
computes the face location of each vertex of a 2D polygon on tmesh . | |
template<class K , class TriangleMesh > | |
std::vector< std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > > | CGAL::Vector_graphics_on_surfaces::trace_geodesic_polygons (const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > ¢er, const std::vector< std::vector< typename K::Point_2 > > &polygons, const typename K::FT scaling, const TriangleMesh &tmesh, const Dual_geodesic_solver< typename K::FT > &solver={}) |
computes for each vertex of each polygon in polygons a face location on tmesh , where center represents the center of the 2D bounding box of the polygons. | |
template<class K , class TriangleMesh > | |
std::vector< std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > > | CGAL::Vector_graphics_on_surfaces::trace_geodesic_label (const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > ¢er, const std::vector< std::vector< typename K::Point_2 > > &polygons, const typename K::FT scaling, const TriangleMesh &tmesh, const Dual_geodesic_solver< typename K::FT > &solver={}) |
computes for each vertex of each polygon in polygons a face location on tmesh , where center represents the center of the 2D bounding box of the polygons. | |
template<class K , class TriangleMesh > | |
std::vector< std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > > | CGAL::Vector_graphics_on_surfaces::trace_geodesic_label_along_curve (const std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > &supporting_curve, const std::vector< std::vector< typename K::Point_2 > > &polygons, const typename K::FT scaling, const typename K::FT padding, const bool is_centered, const TriangleMesh &tmesh, const Dual_geodesic_solver< typename K::FT > &solver={}) |
computes for each vertex of each polygon in polygons a face location on tmesh along the curve supporting_curve . | |
template<class K , class TriangleMesh > | |
std::vector< std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > > | CGAL::Vector_graphics_on_surfaces::trace_bezier_curves (const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > ¢er, const std::vector< std::array< typename K::Vector_2, 4 > > &directions, const std::vector< std::array< typename K::FT, 4 > > &lengths, const int num_subdiv, const TriangleMesh &tmesh, const Dual_geodesic_solver< typename K::FT > &solver={}) |
computes a path representing a Bézier curve defined by four control points. | |
template<class K , class TriangleMesh > | |
std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > | CGAL::Vector_graphics_on_surfaces::trace_polyline_of_bezier_curves (const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > ¢er, const std::vector< typename K::Vector_2 > &directions, const std::vector< typename K::FT > &lengths, bool is_closed, const int num_subdiv, const TriangleMesh &tmesh, const Dual_geodesic_solver< typename K::FT > &solver={}) |
computes a path representing a Bézier polyline (a sequence of Bézier curves having an identical control points, that is the fourth control point of the nth curve is the first control point of the (n+1)th curve). | |
template<class K , class TriangleMesh , class VNM , class FNM , class OutputIterator > | |
OutputIterator | CGAL::Vector_graphics_on_surfaces::refine_mesh_along_paths (const std::vector< std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > > &paths, TriangleMesh &tmesh, VNM vnm, FNM fnm, OutputIterator out) |
refines tmesh so that each path in paths corresponds to a set edges of tmesh after the call. | |
void CGAL::Vector_graphics_on_surfaces::locally_shortest_path | ( | CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, FT > | src, |
CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, FT > | tgt, | ||
const TriangleMesh & | tmesh, | ||
EdgeLocationRange & | edge_locations, | ||
const Dual_geodesic_solver< FT > & | solver | ||
) |
#include <CGAL/Vector_graphics_on_surfaces/locally_shortest_path.h>
computes an approximated geodesic shortest path between two locations on a tmesh
.
The pointssrc
and tgt
must be on the same connected component.
TriangleMesh | a model of FaceListGraph and EdgeListGraph |
FT | floating point number type (float or double) |
EdgeLocationRange | a model of BackInsertionSequence whose value type CGAL::Polygon_mesh_processing::Edge_location<FT> . |
src | source of the path |
tgt | target of the path |
tmesh | input triangle mesh to compute the path on |
edge_locations | contains the path as a sequence of edge locations. Two consecutives edges e_k and e_kp1 stored in edge_locations are such that face(halfedge(e_k, tmesh), tmesh) == face(opposite(halfedge(e_kp1, tmesh), tmesh), tmesh)) . In parcular, it means that if the path goes through a vertex of tmesh , several edge locations will be reported to maintain this property. Additionally, if src is in the interior of a face f , then the first edge location e_0 of edge_locations is such that f == face(opposite(halfedge(e_0, tmesh), tmesh), tmesh)) . Similarly, if tgt is in the interior of a face f , then the last edge location e_n of edge_locations is such that f == face(halfedge(e_n, tmesh), tmesh) . |
solver | container for the precomputed information. If not initialized, it will be initialized internally. |
add named parameters
should we have halfedge location instead?
std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, FT > > CGAL::Vector_graphics_on_surfaces::recursive_de_Casteljau | ( | const TriangleMesh & | mesh, |
const Bezier_segment< TriangleMesh, FT > & | control_points, | ||
const int | num_subdiv, | ||
const Dual_geodesic_solver< FT > & | solver = Dual_geodesic_solver<FT>() |
||
) |
#include <CGAL/Vector_graphics_on_surfaces/locally_shortest_path.h>
computes a discretization of a Bézier segment defined by the location of four control points on tmesh
.
All control points must be on the same connected component. This functions applies several iterations of the de Casteljau algorithm, and geodesic shortest paths are drawn between the control points.
TriangleMesh | a model of FaceListGraph and EdgeListGraph |
FT | floating point number type (float or double) |
mesh | input triangle mesh to compute the path on |
control_points | control points of the Bézier segment |
num_subdiv | the number of iterations of the subdivision algorithm |
solver | container for the precomputed information. If not initialized, it will be initialized internally. |
add named parameters
do we want to also have a way to return Bézier segments? The output is actually Bézier segments subdivided.
OutputIterator CGAL::Vector_graphics_on_surfaces::refine_mesh_along_paths | ( | const std::vector< std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > > & | paths, |
TriangleMesh & | tmesh, | ||
VNM | vnm, | ||
FNM | fnm, | ||
OutputIterator | out | ||
) |
#include <CGAL/Vector_graphics_on_surfaces/locally_shortest_path.h>
refines tmesh
so that each path in paths
corresponds to a set edges of tmesh
after the call.
Note that each path must be such that for two consecutive face locations, there exists a face in tmesh
containing the two corresponding points.
TriangleMesh | a model of MutableFaceGraph |
K | a model of Kernel with K::FT being a floating point number type (float or double) |
VNM | a model of ReadWritePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Kernel::Vector_3 as value type. |
FNM | a model of ReadWritePropertyMap with boost::graph_traits<TriangleMesh>::face_descriptor as key type and Kernel::Vector_3 as value type. |
OutputIterator | an output iterator accepting boost::graph_traits<TriangleMesh>::halfedge_descriptor to be put |
tmesh | the triangle mesh to be refined |
paths | a path described as a range of edge locations, with the property that for two consecutive edge locations, there exists a face in tmesh containing the two corresponding points. |
vnm | property map associating a normal to each vertex of tmesh that is updated by this function |
fnm | property map associating a normal to each face of tmesh that is updated by this function |
out | output iterator where created halfedges are put |
add named parameters
vnm and fnm are optional
intersection between path or self-intersections are not handle. Should it be? If not what do we do?
out should contain edges, and also existing edges already part of a path
shall we also have the edges in the order of the input rather than all at once
std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > CGAL::Vector_graphics_on_surfaces::straightest_geodesic | ( | const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > & | src, |
const typename K::Vector_2 & | dir, | ||
const typename K::FT | len, | ||
const TriangleMesh & | tmesh | ||
) |
#include <CGAL/Vector_graphics_on_surfaces/locally_shortest_path.h>
computes a path on a triangle mesh that is computed by starting a walk on tmesh
given a direction and a maximum distance.
The distance will not be achieved if a border edge is reached before.
TriangleMesh | a model of FaceListGraph and EdgeListGraph |
FT | floating point number type (float or double) |
tmesh | input triangle mesh to compute the path on |
src | the source of the path |
len | the distance to walk along the straightest |
dir | the initial direction of the walk, given as a 2D vector in the face of src, the halfedge of the face being the y-axis. |
src
) add named parameters
do we want to also have a way to return Bézier segments? The output is actually Bézier segments subdivided.
offer something better than a 2D vector for the direction
std::vector< std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > > CGAL::Vector_graphics_on_surfaces::trace_bezier_curves | ( | const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > & | center, |
const std::vector< std::array< typename K::Vector_2, 4 > > & | directions, | ||
const std::vector< std::array< typename K::FT, 4 > > & | lengths, | ||
const int | num_subdiv, | ||
const TriangleMesh & | tmesh, | ||
const Dual_geodesic_solver< typename K::FT > & | solver = {} |
||
) |
#include <CGAL/Vector_graphics_on_surfaces/locally_shortest_path.h>
computes a path representing a Bézier curve defined by four control points.
Control points are defined by the endpoints of straightest geodesic curves starting from center
along given directions and distances. The iterative de Casteljau subdivision algorithm is applied to create more control points that are then connected with locally shortest paths. The output path is such that for two consecutive face locations, there must exist a face in tmesh
containing the two corresponding points.
TriangleMesh | a model of FaceListGraph and EdgeListGraph |
K | a model of Kernel with K::FT being a floating point number type (float or double) |
center | the location on tmesh where straightest geodesic for the placement of control points starts. The y-axis used is halfedge(center.first, tmesh) . |
directions | contains the direction of the straightest geodesic for each control point |
lengths | contains the length of the straightest geodesic for each control point |
num_subdiv | the number of iterations of the de Casteljau subdivision algorithm |
tmesh | input triangle mesh supporting the vertices of the output polygon |
solver | container for the precomputed information. If not initialized, it will be initialized internally. |
add named parameters
polygon orientation is not handled in the function and should be done outside of the function for now
for better rendering we can group polygons to have one center for the same group of polygon (useful for letters that are not simply connected)
what if boundary is reached
std::vector< std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > > CGAL::Vector_graphics_on_surfaces::trace_geodesic_label | ( | const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > & | center, |
const std::vector< std::vector< typename K::Point_2 > > & | polygons, | ||
const typename K::FT | scaling, | ||
const TriangleMesh & | tmesh, | ||
const Dual_geodesic_solver< typename K::FT > & | solver = {} |
||
) |
#include <CGAL/Vector_graphics_on_surfaces/locally_shortest_path.h>
computes for each vertex of each polygon in polygons
a face location on tmesh
, where center
represents the center of the 2D bounding box of the polygons.
This method starts by considering the segment splitting in two halves along the y-axis the bounding box of the polygons. 2D centers for each polygon are computed on this segment as the intersection with the line splitting the bounding box of the polygon in two halves along the x-axis. The splitting segment is then drawn on tmesh
and the face location of the 2D centers is found. trace_geodesic_polygon()
is then called for each polygon and center, with appropriate directions and distances to have a consistent orientation for the polygons.
TriangleMesh | a model of FaceListGraph and EdgeListGraph |
K | a model of Kernel with K::FT being a floating point number type (float or double) |
center | the location on tmesh corresponding to the center of the 2D bounding box of the polygons. |
polygons | 2D polygons |
scaling | a scaling factor to scale the polygons on tmesh (considering geodesic distances on tmesh ) |
tmesh | input triangle mesh supporting the vertices of the output polygon |
solver | container for the precomputed information. If not initialized, it will be initialized internally. |
add named parameters
polygon orientation is not handled in the function and should be done outside of the function for now
for better rendering we can group polygons to have one center for the same group of polygon (useful for letters that are not simply connected)
what if boundary is reached
std::vector< std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > > CGAL::Vector_graphics_on_surfaces::trace_geodesic_label_along_curve | ( | const std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > & | supporting_curve, |
const std::vector< std::vector< typename K::Point_2 > > & | polygons, | ||
const typename K::FT | scaling, | ||
const typename K::FT | padding, | ||
const bool | is_centered, | ||
const TriangleMesh & | tmesh, | ||
const Dual_geodesic_solver< typename K::FT > & | solver = {} |
||
) |
#include <CGAL/Vector_graphics_on_surfaces/locally_shortest_path.h>
computes for each vertex of each polygon in polygons
a face location on tmesh
along the curve supporting_curve
.
This method starts by considering the segment splitting in two halves along the y-axis the bounding box of the polygons. 2D centers for each polygon are computed on this segment as the intersection with the line splitting the bounding box of the polygon in two halves along the x-axis. The splitting segment is then mapped onto supporting_curve
by first scaling it using scaling
, and using padding
and is_centered
. Face locations of the 2D center are found on supporting_curve
. trace_geodesic_polygon()
is then called for each polygon and center, with appropriate directions and distances to have a consistent orientation for the polygons.
TriangleMesh | a model of FaceListGraph and EdgeListGraph |
K | a model of Kernel with K::FT being a floating point number type (float or double) |
supporting_curve | a path on tmesh that will support the center of the bounding box of each polygon. For two consecutive face locations, there must exist a face in tmesh containing the two corresponding points. |
polygons | 2D polygons |
scaling | a scaling factor to scale the polygons on tmesh (considering geodesic distances on tmesh ) |
padding | padding applied at the beggining of supporting curve to start the drawing |
is_centered | is true , padding is ignored and the bounding box of the polygon is centered on the supporting curve (in 1D) |
tmesh | input triangle mesh supporting the vertices of the output polygon |
solver | container for the precomputed information. If not initialized, it will be initialized internally. |
add named parameters
polygon orientation is not handled in the function and should be done outside of the function for now
for better rendering we can group polygons to have one center for the same group of polygon (useful for letters that are not simply connected)
check padding is ignored if is_centered is used + update the doc if not
doc what happens if supporting curve is not long enough + boundary reached
std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > CGAL::Vector_graphics_on_surfaces::trace_geodesic_polygon | ( | const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > & | center, |
const std::vector< typename K::Vector_2 > & | directions, | ||
const std::vector< typename K::FT > & | lengths, | ||
const TriangleMesh & | tmesh, | ||
const Dual_geodesic_solver< typename K::FT > & | solver = {} |
||
) |
#include <CGAL/Vector_graphics_on_surfaces/locally_shortest_path.h>
computes the face location of each vertex of a 2D polygon on tmesh
.
The vertices of the polygon are given as a pair of direction and distance with respect to a point corresponding to center
on tmesh
.
TriangleMesh | a model of FaceListGraph and EdgeListGraph |
K | a model of Kernel with K::FT being a floating point number type (float or double) |
center | the location on tmesh used as reference. The y-axis used for coordinates is halfedge(center.first, tmesh) . |
directions | contains the direction one need to move from center to reach each vertex of the polygon. |
lengths | the distance one need to move from center along the direction at the same position in directions to reach each vertex of the polygon. |
tmesh | input triangle mesh supporting the vertices of the output polygon |
solver | container for the precomputed information. If not initialized, it will be initialized internally. |
add named parameters
polygon orientation is not handled in the function and should be done outside of the function for now
offer something better than a 2D vector for the direction
directly handle polyline?
why the first polygon vertex is duplicated by the function? (most probably for the example but it shouldn't be done here)
std::vector< std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > > CGAL::Vector_graphics_on_surfaces::trace_geodesic_polygons | ( | const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > & | center, |
const std::vector< std::vector< typename K::Point_2 > > & | polygons, | ||
const typename K::FT | scaling, | ||
const TriangleMesh & | tmesh, | ||
const Dual_geodesic_solver< typename K::FT > & | solver = {} |
||
) |
#include <CGAL/Vector_graphics_on_surfaces/locally_shortest_path.h>
computes for each vertex of each polygon in polygons
a face location on tmesh
, where center
represents the center of the 2D bounding box of the polygons.
This method computes the location of the center of the bounding box of each polygon on the mesh with respect to center
and calls trace_geodesic_polygon()
with that center with appropriate directions and distances to have a consistent orientation for the polygons.
TriangleMesh | a model of FaceListGraph and EdgeListGraph |
K | a model of Kernel with K::FT being a floating point number type (float or double) |
center | the location on tmesh corresponding to the center of the 2D bounding box of the polygons. |
polygons | 2D polygons |
scaling | a scaling factor to scale the polygons on tmesh (considering geodesic distances on tmesh ) |
tmesh | input triangle mesh supporting the vertices of the output polygon |
solver | container for the precomputed information. If not initialized, it will be initialized internally. |
add named parameters
polygon orientation is not handled in the function and should be done outside of the function for now
for better rendering we can group polygons to have one center for the same group of polygon (useful for letters that are not simply connected)
what if boundary is reached
std::vector< CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > > CGAL::Vector_graphics_on_surfaces::trace_polyline_of_bezier_curves | ( | const CGAL::Polygon_mesh_processing::Face_location< TriangleMesh, typename K::FT > & | center, |
const std::vector< typename K::Vector_2 > & | directions, | ||
const std::vector< typename K::FT > & | lengths, | ||
bool | is_closed, | ||
const int | num_subdiv, | ||
const TriangleMesh & | tmesh, | ||
const Dual_geodesic_solver< typename K::FT > & | solver = {} |
||
) |
#include <CGAL/Vector_graphics_on_surfaces/locally_shortest_path.h>
computes a path representing a Bézier polyline (a sequence of Bézier curves having an identical control points, that is the fourth control point of the nth curve is the first control point of the (n+1)th curve).
Control points are defined by the endpoints of straightest geodesic curves starting from center
along given directions and distances. The iterative de Casteljau subdivision algorithm is applied to create more control points that are then connected with locally shortest paths. The output path is such that for two consecutive face locations, there must exist a face in tmesh
containing the two corresponding points. The first portion of the curve is defined by the 4 first values in directions
and lengths
. The second portion is defined by the 4'th value and the next 3, and so on until the end is reached. If is_closed
is true, then the first value will be used with the last three to define the last portion.
TriangleMesh | a model of FaceListGraph and EdgeListGraph |
K | a model of Kernel with K::FT being a floating point number type (float or double) |
center | the location on tmesh where straightest geodesic for the placement of control points starts. The y-axis used is halfedge(center.first, tmesh) . |
directions | contains the direction of the straightest geodesic for each control point |
lengths | contains the length of the straightest geodesic for each control point |
num_subdiv | the number of iterations of the de Casteljau subdivision algorithm |
is_closed | if true [directions/lengths].front() will be used as additional last point, generating a closed path |
tmesh | input triangle mesh supporting the vertices of the output polygon |
solver | container for the precomputed information. If not initialized, it will be initialized internally. |
add named parameters
polygon orientation is not handled in the function and should be done outside of the function for now
for better rendering we can group polygons to have one center for the same group of polygon (useful for letters that are not simply connected)
what if boundary is reached