ddg.halfedge.delaunay module
Module delaunay gives methods at hand to check for Delaunay property of
a triangulation and to transform it into a Delaunay triangulation
via edge-flips.
- ddg.halfedge.delaunay.calculate_angle(e)[source]
Calculates the inner angle between an edge and its consecutive next edge within a face.
- Parameters:
- eedge
The edge pointing towards the inner angle that will be calculated.
- Returns:
- angledouble
The angle between
eande.nex.
- Raises:
- ValueError
If the edge is a boundary edge, i.e. does not belong to a face. If the triangle is degenerated or the triangle inequality doesn’t hold.
- ddg.halfedge.delaunay.opp_halfangle_tan(e)[source]
Calculates the tangent of half the angle that lies opposite the given halfedge in the corresponding triangle face.
- Parameters:
- eedge
The edge whose opposite half-angle tangent will be calculated.
- Returns:
- tan_alpha_kij_halfdouble
The tangent of half the angle that lies opposite edge
ein the corresponding triangle.
- Raises:
- ValueError
If the edge is a boundary edge, i.e. does not belong to a face. If the triangle is degenerated or the triangle inequality doesn’t hold.
Notes
The tangents of the half angles of the triangles are computed using the half-angle theorem. For the angle \(\alpha\) opposite
ain a triangle with side lengtha,b,cit is\[\tan \frac{\alpha}{2} = \sqrt{\frac{(a-b+c)(a+b-c)}{(a+b+c)(-a+b+c)}}.\]
- ddg.halfedge.delaunay.opp_angle_cot(e)[source]
Calculates the cotangent of the angle that is lying opposite the given halfedge in the corresponding triangle face.
- Parameters:
- eedge
The edge whose opposite angle co-tangent will be calculated.
- Returns:
- cot_alphadouble
The cotangent of the angle alpha opposite the given edge e.
Notes
The cotangent is calculated using the half-angle tangents .. math:: cot(alpha) = frac{1-tan^2(alpha /2)}{2tan(alpha /2)}.
These cotangents can be used to form the so called cotan weight \(w_{ij} = \cot(\alpha^k_{ij}) + \cot(\alpha^l_{ij})\) for an edge eij with the adjacent triangles
ijkandilj. An edge is Delaunay, iff its cotan weight is non-negative.
- ddg.halfedge.delaunay.euclidean_cotan_weight(e)[source]
Calculates the Eulcidean cotan weight for an edge of a triangulation.
- Parameters:
- eedge
The edge whose cotan weight shall be calculated.
- Returns:
- wijdouble
The cotan weight of the given edge.
Notes
The Euclidean cotan weight of an edge
eijwith adjacent trianglesijkandiljis \(w_{ij} = \cot(\alpha^k_{ij}) + \cot(\alpha^l_{ij})\).For an edge
eijopposite to the boundary in triangleijkwe set \(w_{ij} = \cot(\alpha^k_{ij})\).The cotangents in the cotan weights are computed using the half angle tangents of the angles opposite the edge that is to be checked.
- ddg.halfedge.delaunay.is_delaunay(edge_or_surface, evaluation_fct=<function evaluate_cotan_weight>)[source]
Checks if an edge adjacent to two triangles or a triangulated surface is Delaunay.
- Parameters:
- edge_or_surfaceedge or surface
The edge or surface that is checked for Delaunay property.
- evaluation_fctfunction(edge) -> bool
A function that takes an edge as argument and determines how being delaunay is evaluated. By default the Euclidean cotan weight is used to evaluate the delaunay property.
- Returns:
- bool
- Raises:
- ValueError
If the input is a surface that is not a triangulation. If the input is an edge whose face is not a triangle.
Notes
By default the Delaunay property is checked by using the Euclidean cotan weights \(w_{ij} = \cot(\alpha^k_{ij}) + \cot(\alpha^l_{ij})\) for an edge eij with adjacent triangles
ijkandilj. For boundary edges,is_delaunay(eij)returns true, for an interior edge it returns true if the cotan weight wij is nonnegative, false otherwise.
- ddg.halfedge.delaunay.update_edgelength(e, length_function)[source]
- Updates the edgelength attribute
lenof a given halfedge and its opposite by use of a given function that determines the new length.
- Parameters:
- eedge
The halfedge whose length attribute
lenwill be updated according to the givenlength_function. The corresponding opposite edgelengthe.opp.lenwill also be updated.- length_functionfunction(edge) -> float
A function that calculates a new edgelength for a given edge.
- Returns:
- None
- Updates the edgelength attribute
- ddg.halfedge.delaunay.intrinsic_lengthfunction(e)[source]
- Function that calculates the new intrinsic edgelength for an edge that will
be flipped.
For calculation of the new edgelength after a flip the distance between the new head and tail coordinates is measured.
- Parameters:
- eedge
The edge whose edgelength after a flip will be calculated in extrinsic manner.
- Returns:
- new_edgelengthfloat
The new intrinsic length of the corresponding flipped edge.
Notes
The new edgelength \(\ell_{lk}\) of the flipped edge is computed using half-angle tangents and law of cosine:
\[\tan(\frac{\alpha + \delta}{2}) = \frac{\tan(\alpha /2) + \tan(delta /2)} {1-\tan(\alpha /2)\tan(\delta /2)}\]then
\[\cos(\alpha + \delta) = \frac{1-\tan^2(\frac{\alpha + \delta}{2})} {1+\tan^2(\frac{\alpha + \delta}{2})}\]and for the new length
\[\ell_{lk} = \sqrt{b^2 + c^2 -2bc\cos(\alpha + \delta)}.\]
- ddg.halfedge.delaunay.extrinsic_lengthfunction(e)[source]
- Function that calculates the new extrinsic edgelength for an edge that will
be flipped.
For calculation of the new edgelength after a flip the distance between the new head and tail coordinates is measured.
- Parameters:
- eedge
The edge whose edgelength after a flip will be calculated in extrinsic manner.
- Returns:
- new_edgelengthfloat
The new extrinsic length of the corresponding flipped edge.
- ddg.halfedge.delaunay.flip(e, length_function=<function intrinsic_lengthfunction>)[source]
Flips an edge within a quadrilateral of two adjacent triangles to the other diagonal edge and by default updates its new edgelength in intrinsic manner.
- Parameters:
- eedge
The edge within a kite of two adjacent triangles that shall be flipped.
- length_functionfunction(edge) -> float
A function that takes an edge (before it is being flipped) as argument and calculates a new edgelength which will be assigned to the flipped edge. By default the Euclidean intrinsic length is used to update the length of
e.
- Returns:
- None
- Raises:
- ValueError
If the edge is not lying between two adjacent triangles.
Notes
For an edge
e_ijbetween the trianglesijkandiljflip(e_ij)will flip both halfedgese_ijande_jito the corresponding other diagonal halfedgese_lkande_klin the quadrilateraliljk. And all references are updated accordingly.This method performs an intrinsic edge flip by default, the flipped edge is a geodesic segment in the corresponding surface thus the shape of the the original embedded mesh does not change. The new edgelength \(\ell_{lk}\) of the flipped edge is computed using half-angle tangents and law of cosine:
\[\tan(\frac{\alpha + \delta}{2}) = \frac{\tan(\alpha /2) + \tan(delta /2)} {1-\tan(\alpha /2)\tan(\delta /2)}\]then
\[\cos(\alpha + \delta) = \frac{1-\tan^2(\frac{\alpha + \delta}{2})} {1+\tan^2(\frac{\alpha + \delta}{2})}\]and for the new length
\[\ell_{lk} = \sqrt{b^2 + c^2 -2bc\cos(\alpha + \delta)}.\]
- ddg.halfedge.delaunay.delaunay_flip_algorithm(original_surface, return_overlay=False, length_function=<function intrinsic_lengthfunction>, evaluation_fct=<function evaluate_cotan_weight>)[source]
Transforms a triangulation into a Delaunay triangulation by flipping edges.
An intrinsic Delaunay edge flipping algorithm that intrinsically flips edges of an initial triangulation such that in the end all edges are Delaunay and so the surface is. One can pass any desired evaluation function for determining the delaunay property as well as a length function for calculation of new edgelength. The result of the intrinsic_delaunay algorithm is called the corresponding overlay structure. This surface may be non-regular, nevertheless we call it a triangulation and do not require strong regularity for our purpose.
- Parameters:
- original_surfacesurface
The triangulated surface that shall be transferred into a Delaunay triangulation via edge flips.
- return_overlaybool
If True, an overlay structure will be created and returned, which stores the intrinsic Delaunay Triangulation with edge attributes
crsgiving the crossings, i.e. a list of edge tuples oforiginal_surface. If False, the original surface will be changed to be Delaunay and returned. False by default.- length_functionfunction(edge) -> float
A function that takes an edge (before it is being flipped) as argument and calculates a new edgelength which will be assigned to the flipped edge. By default the Euclidean intrinsic_lengthfunction is used to update the length of
e.- evaluation_fctfunction(edge) -> bool
A function that takes an edge as argument and determines how being delaunay is evaluated. By default the Euclidean cotan weight is used to evaluate the delaunay property.
- Returns:
- surfsurface
The resulting Delaunay triangulation after edge flipping. Depending on the argument
return_overlayeither the original surface will be transformed and returned or an overlay structure will be created and returned, which stores the intrinsic Delaunay Triangulation with edge attributescrsgiving the crossings.
- Raises:
- ValueError
If the input is a surface that is not a triangulation.
Notes
By default the Delaunay property is evaluated using the corresponding Euclidean cotan weights. If an edge is not Delaunay it will be flipped and its length will be updated, by default using the intrinsic new length.
The usual definition of triangulation, which implies strong regularity, is to narrow for our purpose. Even if the input triangulation is strongly regular, the result of intrinsic edge flipping may be non-regular. Therefore we do not require triangulations to be regular.
In case of a flip, the four neighbouring edges will be stacked for later check.
- ddg.halfedge.delaunay.update_crossings(e)[source]
Updates the crossings of an edge belonging to an overlay structure after its flip, i.e. stores all original edges that are now crossed by this current one.
- Parameters:
- eedge
The edge whose crossings attribute
e.crswill be updated.
- Returns:
- None
Notes
The crossings of an edge, meaning which edges of the original surface are crossed by it, are stored in the edge attribute called
crsas a list of edge-tuples of the form(edge, edge.opp). This list gives the halfedges in order of their crossing frome.tailtoe.head.
- ddg.halfedge.delaunay.init_crossings_attr(surf)[source]
Adds a list as crossings attribute ‘crs’ to all edges of the surface.
In this list all crossings of an edge of the current triangulation with edges of the original triangulation are saved. This is used for bookkeeping of and maintaining the overlay structure during flipping.
- Parameters:
- surfsurface
The surface whose edges will be added a crossings attribute that can be addressed via ‘edge.crs’.
- Returns:
- None
- ddg.halfedge.delaunay.init_overlay_structure(surf)[source]
Creates an inital overlay structure for the given surface.
A copy of the given surface is made and provided with the necessary attributes in order to function as an initial overlay structure.
- Parameters:
- surfsurface
- Returns:
- overlaysurface
The copy of the given surface with reference to the original vertices via the attribute
original_vertexand the edge attributecrsfor storing the crossings of original edges.
- class ddg.halfedge.delaunay.SurfaceWithOverlay(original_surf, overlay_surf)[source]
Bases:
objectA class for constructing a visualisable intrinsic delaunay triangulation from a halfedge surface and its intrinsic delaunay overlay structure.
The overlay structure is the result of performing the intrinsic edge flip algorithm on the initial surface. We call this surface a triangulation, but it may be non-regular. Also, as the word intrinsic implies, its edges are geodesics on the initial surface, so not straight lines in \(\mathbb{R}^3\) anymore. For its visualization we need to take further steps. They are implemented by methods of this class.
- get_intrinsic_delaunay_surf()[source]
Constructs a visualisable intrinsic Delaunay triangulation from a SurfaceWithOverlay.
- Returns:
- delaunay_surfsurf
A visualisable delaunay triangulation corresponding to a SurfaceWithOverlay.
- layout_of_crossed_faces(e)[source]
Finds 2D coordinates for all vertices in a strip of triangles (ORI) that are crossed by the given edge (OLY).
- Parameters:
- eedge
The edge for which the strip of crossed triangles (ORI) will be layed out in the plane.
- Returns:
- uv_mapV(ORI) -> [float, float]
A map assigning each vertex of the original surface a 2D coordinate.
- calculate_crs_coordinates(e, new_verts)[source]
- Calculates and assigns 3D coordinates to the new vertices of a delaunay
surface that correspond to a given edge in the overlay surface.
When we evolve a delaunay surf from an overlay surf, we subdivide those edges that crossed original ones in order to achieve straight segments as edges. Thereby we added new vertices to the delaunay surf which we now calculate coordinates for. First we lay out the strip of original triangles the overlay edge runs through in the plane (2D). Then the coordinates of the new verts (DEL) are calculated as intersections of this overlay edge with the crossed ones (ORI). The intersection will be given by its barycentric coefficients relative to the endpoints of the intersected segments. By this we do not only get the coordinates for the intersection in the plane but immediately also those for the corresponding point on the surface.
- Parameters:
- eedge (OLY - in an overlay surface)
An edge of the overlay surf for which the coordinates of its corresponding crossings vertices shall be calculated for.
- new_vertslist of verts (DEL)
The vertices in the delaunay surf that were added via subdivision of
ewhich now shall be given 3D coordinates.
- Returns:
- None
- ddg.halfedge.delaunay.add_parser(sub_parsers)[source]
Parses arguments for the Delaunay edge flip algo command line application.
Takes the user-desired input arguments of the command line and performs the Delaunay edge flip algo on the given input triangulation with it.
- Parameters:
- sub_parsersArgumentParser
The parent argument subparser of
ddg.main()which this parser shall belong to.
- Returns:
- None
Notes
Performs the algorithm by calling
execute_flip_algo(args)with the parsed arguments.
- ddg.halfedge.delaunay.execute_flip_algo(args)[source]
Performs the Delaunay edge flip algorithm on a given triangulation.
Takes a triangulation given in obj or json format as input and performs the Delaunay edge flip algorithm according to chosen optional arguments. The resulting flipped structure is returned in a new file of format obj (default) or json.
- Parameters:
- argsan object holding all parsed arguments as attributes
The input arguments handed in by a user in the command line.
- args.fileobj or json file
The input triangulation that shall be transformed via the edge flip algorithm, given in obj or json format.
- Returns:
- Writes the resulting structure to a new obj/json file.
Notes
For performing the edge flip algorithm the given input will be converted into halfedge datastructure. And converted back into obj or json format after performance.