CGAL 5.3 - Triangulated Surface Mesh Deformation
CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM > Class Template Reference

#include <CGAL/Surface_mesh_deformation.h>

Definition

template<class TM, class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
class CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >

Class providing the functionalities for deforming a triangulated surface mesh.

Template Parameters
TMa model of HalfedgeGraph
VIMa model of ReadablePropertyMap with Surface_mesh_deformation::vertex_descriptor as key and unsigned int as value type. The default is boost::property_map<TM, boost::vertex_index_t>::type.
HIMa model of ReadablePropertyMap with Surface_mesh_deformation::halfedge_descriptor as key and unsigned int as value type. The default is boost::property_map<TM, boost::halfedge_index_t>::type.
TAGtag for selecting the deformation algorithm
WCa model of SurfaceMeshDeformationWeights, with WC::Triangle_mesh being TM. If TAG is ORIGINAL_ARAP, the weights must be positive to guarantee a correct energy minimization. The default is the cotangent weighting scheme. In case TAG is ORIGINAL_ARAP, negative weights are clamped to zero.
STa model of SparseLinearAlgebraWithFactorTraits_d. If Eigen 3.2 (or greater) is available and CGAL_EIGEN3_ENABLED is defined, then an overload of Eigen_solver_traits is provided as default parameter.
Eigen::SparseLU<
Eigen::COLAMDOrdering<int> > >
CRa model of DeformationClosestRotationTraits_3. If Eigen 3.1 (or greater) is available and CGAL_EIGEN3_ENABLED is defined, Deformation_Eigen_polar_closest_rotation_traits_3 is provided as default parameter.
VPMa model of ReadWritePropertyMap with Surface_mesh_deformation::vertex_descriptor as key and a point as value type. The point type must be a model of RawPoint_3. The default is boost::property_map<TM, CGAL::vertex_point_t>::type.
Examples:
Surface_mesh_deformation/all_roi_assign_example.cpp, Surface_mesh_deformation/custom_weight_for_edges_example.cpp, Surface_mesh_deformation/deform_mesh_for_botsch08_format_sre_arap.cpp, Surface_mesh_deformation/deform_polyhedron_with_custom_pmap_example.cpp, and Surface_mesh_deformation/k_ring_roi_translate_rotate_example.cpp.

Types

typedef TM Triangle_mesh
 Triangle mesh type.
 
typedef VIM Vertex_index_map
 vertex index map type
 
typedef HIM Hedge_index_map
 halfedge index map type
 
typedef WC Weight_calculator
 weight calculator functor type
 
typedef ST Sparse_linear_solver
 sparse linear solver type
 
typedef CR Closest_rotation_traits
 closest rotation traits type
 
typedef VPM Vertex_point_map
 vertex point map type
 
typedef boost::graph_traits< Triangle_mesh >::vertex_descriptor vertex_descriptor
 The type for vertex descriptor.
 
typedef boost::graph_traits< Triangle_mesh >::halfedge_descriptor halfedge_descriptor
 The type for halfedge descriptor.
 
typedef boost::property_traits< Vertex_point_map >::value_type Point
 The 3D point type, model of RawPoint_3
 
typedef std::vector< vertex_descriptorRoi_vertex_range
 A constant iterator range over the vertices of the region-of-interest. More...
 

Construction

 Surface_mesh_deformation (Triangle_mesh &triangle_mesh, Vertex_index_map vertex_index_map=unspecified_internal_vertex_index_map, Hedge_index_map hedge_index_map=unspecified_internal_halfedge_index_map, Vertex_point_map vertex_point_map=get(boost::vertex_point, triangle_mesh), Weight_calculator weight_calculator=Weight_calculator())
 The constructor of a deformation object. More...
 

Preprocessing

void clear_roi_vertices ()
 Erases all the vertices from the region-of-interest (control vertices included).
 
void clear_control_vertices ()
 Erases all the vertices from the set of control vertices.
 
bool insert_control_vertex (vertex_descriptor vd)
 Inserts a vertex in the set of control vertices. More...
 
template<class InputIterator >
void insert_control_vertices (InputIterator begin, InputIterator end)
 Inserts a range of vertices in the set of control vertices. More...
 
bool erase_control_vertex (vertex_descriptor vd)
 Erases a vertex from the set of control vertices. More...
 
template<class InputIterator >
void insert_roi_vertices (InputIterator begin, InputIterator end)
 Inserts a range of vertices in the region-of-interest. More...
 
bool insert_roi_vertex (vertex_descriptor vd)
 Inserts a vertex in the region-of-interest. More...
 
bool erase_roi_vertex (vertex_descriptor vd)
 Erases a vertex from the region-of-interest and the set of control vertices. More...
 
bool preprocess ()
 Preprocessing function that need to be called each time the region-of-interest or the set of control vertices are changed before calling deform(). More...
 

Deformation

void set_target_position (vertex_descriptor vd, const Point &target_position)
 Sets the target position of a control vertex. More...
 
template<class Vect >
void translate (vertex_descriptor vd, const Vect &t)
 Updates the target position of vd by applying a translation of vector t. More...
 
template<class InputIterator , class Vect >
void translate (InputIterator begin, InputIterator end, const Vect &t)
 Equivalent to calling the overload taking only one control vertex, for each vertex in the range [begin,end[. More...
 
template<typename Quaternion , typename Vect >
void rotate (vertex_descriptor vd, const Vect &to_rotation_center, const Quaternion &quat)
 Updates the target position of vd by applying to its last target position a rotation defined by the quaternion quat, the center of the rotation being the origin translated by to_rotation_center . More...
 
template<typename InputIterator , typename Vect , typename Quaternion >
void rotate (InputIterator begin, InputIterator end, const Vect &to_rotation_center, const Quaternion &quat)
 Equivalent to calling the overload taking only one control vertex, for each vertex in the range [begin,end[. More...
 
const Pointtarget_position (vertex_descriptor vd)
 Returns the target position of a control vertex. More...
 
void deform ()
 Deforms the region-of-interest according to the deformation algorithm, using the target positions of each control vertex set by using rotate(), translate(), or set_target_position(). More...
 
void deform (unsigned int iterations, double tolerance)
 Same as deform() but the number of iterations and the tolerance are one-time parameters. More...
 
void reset ()
 Resets the points associated to the vertices of the region-of-interest at their initial positions at time of the functor construction or after the last call to overwrite_initial_geometry(). More...
 
void overwrite_initial_geometry ()
 Sets the initial positions of the vertices from the region-of-interest to the current positions. More...
 

Utilities

unsigned int iterations ()
 Gets the default number of iterations (5) or the value passed to the function set_iterations()
 
double tolerance ()
 Gets the default tolerance parameter (1e-4) or the value passed to the function set_tolerance()
 
void set_iterations (unsigned int iterations)
 Sets the maximum number of iterations ran by deform()
 
void set_tolerance (double tolerance)
 Sets the tolerance of the convergence used in deform(). More...
 
const Roi_vertex_rangeroi_vertices () const
 Returns the range of vertices in the region-of-interest.
 
bool is_roi_vertex (vertex_descriptor vd) const
 Tests whether a vertex is inside the region-of-interest. More...
 
bool is_control_vertex (vertex_descriptor vd) const
 Tests whether a vertex is a control vertex. More...
 
const Triangle_meshtriangle_mesh () const
 Provides access to the triangle mesh being deformed.
 
void set_sre_arap_alpha (double a)
 Sets the alpha coefficient that determines the weight of the bending term (rotation smoothness) for the SRE-ARAP deformation technique. More...
 

Member Typedef Documentation

◆ Roi_vertex_range

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
typedef std::vector<vertex_descriptor> CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::Roi_vertex_range

A constant iterator range over the vertices of the region-of-interest.

It is a model of ConstRange with vertex_descriptor as iterator value type.

Constructor & Destructor Documentation

◆ Surface_mesh_deformation()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::Surface_mesh_deformation ( Triangle_mesh triangle_mesh,
Vertex_index_map  vertex_index_map = unspecified_internal_vertex_index_map,
Hedge_index_map  hedge_index_map = unspecified_internal_halfedge_index_map,
Vertex_point_map  vertex_point_map = get(boost::vertex_point, triangle_mesh),
Weight_calculator  weight_calculator = Weight_calculator() 
)

The constructor of a deformation object.

Precondition
triangle_mesh consists of only triangular facets
Parameters
triangle_meshtriangulated surface mesh to deform
vertex_index_mapa property map which associates a unique id to each vertex, between 0 to num_vertices(triangle_mesh)-1.
hedge_index_mapproperty map which associates a unique id to each halfedge, between 0 to 2*num_edges(triangle_mesh)-1.
vertex_point_mapproperty map which associates a point to each vertex of the triangle mesh.
weight_calculatorfunction object or pointer for weight calculation

Member Function Documentation

◆ deform() [1/2]

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::deform ( )

Deforms the region-of-interest according to the deformation algorithm, using the target positions of each control vertex set by using rotate(), translate(), or set_target_position().

The points associated to each vertex of the input triangle mesh that are inside the region-of-interest are updated.

Note
Nothing happens if preprocess() returns false.
See also
set_iterations(unsigned int iterations), set_tolerance(double tolerance), deform(unsigned int iterations, double tolerance)

◆ deform() [2/2]

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::deform ( unsigned int  iterations,
double  tolerance 
)

Same as deform() but the number of iterations and the tolerance are one-time parameters.

Parameters
iterationsnumber of iterations for optimization procedure
tolerancetolerance of convergence (see explanations set_tolerance(double tolerance))

◆ erase_control_vertex()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
bool CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::erase_control_vertex ( vertex_descriptor  vd)

Erases a vertex from the set of control vertices.

Parameters
vdthe vertex to be erased
Returns
true if vd was a control vertex.

◆ erase_roi_vertex()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
bool CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::erase_roi_vertex ( vertex_descriptor  vd)

Erases a vertex from the region-of-interest and the set of control vertices.

Note
At the next call to preprocess(), any vertex that is no longer in the region-of-interest will be assigned to its original position (that is the position of the vertex at the time of construction or after the last call to overwrite_initial_geometry()).
Parameters
vdthe vertex to be erased
Returns
true vd was a vertex from the region-of-interest.

◆ insert_control_vertex()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
bool CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::insert_control_vertex ( vertex_descriptor  vd)

Inserts a vertex in the set of control vertices.

The vertex is also inserted in the region-of-interest if it is not already in it.

Parameters
vdthe vertex to be inserted
Returns
true if vd is not already a control vertex.

◆ insert_control_vertices()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
template<class InputIterator >
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::insert_control_vertices ( InputIterator  begin,
InputIterator  end 
)

Inserts a range of vertices in the set of control vertices.

The vertices are also inserted in the region-of-interest if they are not already in it.

Template Parameters
InputIteratorinput iterator type with vertex_descriptor as value type
Parameters
beginfirst iterator of the range of vertices
endpast-the-end iterator of the range of vertices

◆ insert_roi_vertex()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
bool CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::insert_roi_vertex ( vertex_descriptor  vd)

Inserts a vertex in the region-of-interest.

Parameters
vdthe vertex to be inserted
Returns
true if vd is not already in the region-of-interest.

◆ insert_roi_vertices()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
template<class InputIterator >
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::insert_roi_vertices ( InputIterator  begin,
InputIterator  end 
)

Inserts a range of vertices in the region-of-interest.

Template Parameters
InputIteratorinput iterator with vertex_descriptor as value type
Parameters
beginfirst iterator of the range of vertices
endpast-the-end iterator of the range of vertices

◆ is_control_vertex()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
bool CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::is_control_vertex ( vertex_descriptor  vd) const

Tests whether a vertex is a control vertex.

Parameters
vdthe query vertex
Returns
true if vd has been inserted to (and not erased from) the set of control vertices.

◆ is_roi_vertex()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
bool CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::is_roi_vertex ( vertex_descriptor  vd) const

Tests whether a vertex is inside the region-of-interest.

Parameters
vdthe query vertex
Returns
true if vd has been inserted to (and not erased from) the region-of-interest.

◆ overwrite_initial_geometry()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::overwrite_initial_geometry ( )

Sets the initial positions of the vertices from the region-of-interest to the current positions.

Calling this function has the same effect as creating a new deformation object with the current deformed triangle mesh, keeping the region-of-interest and the set of control vertices.

Note
if the region-of-interest or the set of control vertices have been modified since the last call to preprocess(), it will be called prior to the overwrite.

Advanced

This function might have a non-negligible effect on the result. The Laplacian matrix of the free vertices and the optimal rotations are computed using the original positions of the points associated to the vertices. Thus, if a deformed version of the surface mesh is used as reference, the surface mesh properties the algorithm tries to preserve are those of an altered version (which are already degraded).

◆ preprocess()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
bool CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::preprocess ( )

Preprocessing function that need to be called each time the region-of-interest or the set of control vertices are changed before calling deform().

If not already done, deform() first calls this function.

Advanced

Collects the vertices not in the region-of-interest that are adjacent to a vertex from the region-of-interest (these vertices are internally considered as fixed control vertices). Then assembles and factorizes the Laplacian matrix used in the function deform().

Note
A modification of the set of control vertices or the region-of-interest invalidates the preprocessing data.
Returns
true if successful. A common reason for failure is that the system is rank deficient, which happens for example when all the vertices are in the region-of-interest and no control vertices are set, or if the weighting scheme used features too many zero and breaks the connectivity information.

◆ reset()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::reset ( )

Resets the points associated to the vertices of the region-of-interest at their initial positions at time of the functor construction or after the last call to overwrite_initial_geometry().

Note
if the region-of-interest or the set of control vertices have been modified since the last call to preprocess(), it will be called prior to the reset.

◆ rotate() [1/2]

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
template<typename Quaternion , typename Vect >
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::rotate ( vertex_descriptor  vd,
const Vect &  to_rotation_center,
const Quaternion &  quat 
)

Updates the target position of vd by applying to its last target position a rotation defined by the quaternion quat, the center of the rotation being the origin translated by to_rotation_center .

Template Parameters
Quaternionis a quaternion class with Vect operator*(Quaternion, Vect) returning the product of a quaternion with a vector
Vectis a 3D vector class, with Vect(double x,double y, double z) being a constructor from its Cartesian coordinates and double Vect::operator[](int i) with i=0,1 or 2 returning its Cartesian coordinates.
Parameters
vda control vertex
to_rotation_centerthe vector to translate the origin to the center of the rotation
quatquaternion of the rotation
Precondition
is_control_vertex(vd)
quad represents a rotation

◆ rotate() [2/2]

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
template<typename InputIterator , typename Vect , typename Quaternion >
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::rotate ( InputIterator  begin,
InputIterator  end,
const Vect &  to_rotation_center,
const Quaternion &  quat 
)

Equivalent to calling the overload taking only one control vertex, for each vertex in the range [begin,end[.

Template Parameters
InputIteratorinput iterator type with vertex_descriptor as value type
Quaternionis a quaternion class with Vect operator*(Quaternion, Vect) returning the product of a quaternion with a vector
Vectis a 3D vector class, with Vect(double x,double y, double z) being a constructor from its Cartesian coordinates and double Vect::operator[](int i) with i=0,1 or 2 returning its Cartesian coordinates.
Parameters
beginfirst iterator of the range of vertices
endpast-the-end iterator of the range of vertices
to_rotation_centerthe vector to translate the origin to the center of the rotation
quatquaternion of the rotation
Precondition
quad represents a rotation

◆ set_sre_arap_alpha()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::set_sre_arap_alpha ( double  a)

Sets the alpha coefficient that determines the weight of the bending term (rotation smoothness) for the SRE-ARAP deformation technique.

The range of values can be from 0 to infinity. When alpha=0, the method reverts to ARAP. When alpha is increased, neighboring rotations become similar to each other, where alpha=infinity results in a global rigid transformation of the ROI. The value of alpha depends on the surface area and shape. Since alpha is not too sensitive, it can be tweaked in most cases by multiplying it by powers of 10. The default value for alpha is 0.02.

◆ set_target_position()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::set_target_position ( vertex_descriptor  vd,
const Point target_position 
)

Sets the target position of a control vertex.

Parameters
vdthe control vertex the target position is set
target_positionthe new target position

◆ set_tolerance()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::set_tolerance ( double  tolerance)

Sets the tolerance of the convergence used in deform().

If tolerance==0, no energy based termination criteria is used (preventing to do the energy computation at each iteration step)

tolerance > \(|\mathrm{energy}(m_i) - \mathrm{energy}(m_{i-1})| / \mathrm{energy}(m_i)\) will be used as a termination criterium.

◆ target_position()

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
const Point& CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::target_position ( vertex_descriptor  vd)

Returns the target position of a control vertex.

Parameters
vda control vertex
Precondition
is_control_vertex(vd)

◆ translate() [1/2]

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
template<class Vect >
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::translate ( vertex_descriptor  vd,
const Vect &  t 
)

Updates the target position of vd by applying a translation of vector t.

Template Parameters
Vectis a 3D vector class, with Vect(double x,double y, double z) being a constructor from its Cartesian coordinates and double Vect::operator[](int i) with i=0,1 or 2 returning its Cartesian coordinates.
Parameters
vda control vertex
ttranslation vector
Precondition
is_control_vertex(vd)

◆ translate() [2/2]

template<class TM , class VIM = Default, class HIM = Default, Deformation_algorithm_tag TAG = SPOKES_AND_RIMS, class WC = Default, class ST = Default, class CR = Default, class VPM = Default>
template<class InputIterator , class Vect >
void CGAL::Surface_mesh_deformation< TM, VIM, HIM, TAG, WC, ST, CR, VPM >::translate ( InputIterator  begin,
InputIterator  end,
const Vect &  t 
)

Equivalent to calling the overload taking only one control vertex, for each vertex in the range [begin,end[.

Template Parameters
InputIteratorinput iterator type with vertex_descriptor as value type
Vectis a 3D vector class, with Vect(double x,double y, double z) being a constructor from its Cartesian coordinates and double Vect::operator[](int i) with i=0,1 or 2 returning its Cartesian coordinates.
Parameters
beginfirst iterator of the range of vertices
endpast-the-end iterator of the range of vertices
ttranslation vector