Skip to content
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
3 changes: 3 additions & 0 deletions go.mod
Original file line number Diff line number Diff line change
@@ -0,0 +1,3 @@
module github.com/pradeep-pyro/triangle

go 1.17
5 changes: 5 additions & 0 deletions triangle.go
Original file line number Diff line number Diff line change
Expand Up @@ -3,6 +3,7 @@ package triangle
import (
"fmt"
"strconv"
"unsafe"
)

type SegmentSplitting uint8
Expand Down Expand Up @@ -116,6 +117,8 @@ func ConstrainedDelaunay(pts [][2]float64, segs [][2]int32,
in := NewTriangulateIO()
out := NewTriangulateIO()
defer FreeTriangulateIO(in)
// this is duplicated in out, so free it only once
defer trifree(unsafe.Pointer(in.ct.holelist))
defer FreeTriangulateIO(out)

in.SetPoints(pts)
Expand All @@ -139,6 +142,8 @@ func ConformingDelaunay(pts [][2]float64, segs [][2]int32,
in := NewTriangulateIO()
out := NewTriangulateIO()
defer FreeTriangulateIO(in)
// this is duplicated in out, so free it only once
defer trifree(unsafe.Pointer(in.ct.holelist))
defer FreeTriangulateIO(out)

in.SetPoints(pts)
Expand Down
93 changes: 93 additions & 0 deletions triangle_test.go
Original file line number Diff line number Diff line change
@@ -0,0 +1,93 @@
package triangle

import "testing"

// To test if this library adheres to the passing pointers rules as mentioned in
// https://pkg.go.dev/cmd/cgo#hdr-Passing_pointers, run the test suite with the
// command
//
// GODEBUG=cgocheck=2 go test -gcflags=all=-d=checkptr

// Points forming a spiral shape (https://www.cs.cmu.edu/~quake/spiral.node)
var delaunayVoronoiTestPts = [][2]float64{
{0, 0}, {-0.416, 0.909}, {-1.35, 0.436}, {-1.64, 0.549},
{-1.31, -1.51}, {-0.532, -2.17}, {0.454, -2.41}, {1.45, -2.21},
{2.29, -1.66}, {2.88, -0.838}, {3.16, 0.131}, {3.12, 1.14},
{2.77, 2.08}, {2.16, 2.89}, {1.36, 3.49},
}

func TestDelaunay(t *testing.T) {
triangles := Delaunay(delaunayVoronoiTestPts)
const numExpectedTriangles = 16
assertEq(t, numExpectedTriangles, len(triangles),
"Did not get expected number of triangles")
}

func TestVoronoi(t *testing.T) {
vertices, edges, rayOrigins, rayDirections := Voronoi(delaunayVoronoiTestPts)

const (
numExpectedVertices = 16
numExpectedEdges = numExpectedVertices + 2
numExpectedRayOrigins = 12
numExpectedRayDirections = numExpectedRayOrigins
)

assertEq(t, numExpectedVertices, len(vertices),
"Did not get expected number of vertices")
assertEq(t, numExpectedEdges, len(edges),
"Did not get expected number of edges")
assertEq(t, numExpectedRayOrigins, len(rayOrigins),
"Did not get expected number of ray origins")
assertEq(t, numExpectedRayDirections, len(rayDirections),
"Did not get expected number of ray directions")
}

func TestConstrainedConformingDelaunay(t *testing.T) {
// Points forming the shape of letter "A"
var pts = [][2]float64{{0.200000, -0.776400}, {0.220000, -0.773200},
{0.245600, -0.756400}, {0.277600, -0.702000}, {0.488800, -0.207600}, {0.504800, -0.207600}, {0.740800, -0.7396}, {0.756000, -0.761200},
{0.774400, -0.7724}, {0.800000, -0.776400}, {0.800000, -0.792400}, {0.579200, -0.792400}, {0.579200, -0.776400}, {0.621600, -0.771600},
{0.633600, -0.762800}, {0.639200, -0.744400}, {0.620800, -0.684400}, {0.587200, -0.604400}, {0.360800, -0.604400}, {0.319200, -0.706800},
{0.312000, -0.739600}, {0.318400, -0.761200}, {0.334400, -0.771600}, {0.371200, -0.776400}, {0.371200, -0.792400}, {0.374400, -0.570000},
{0.574400, -0.5700}, {0.473600, -0.330800}, {0.200000, -0.792400},
}
// Segments connecting the points
var segs = [][2]int32{{28, 0}, {0, 1}, {1, 2}, {2, 3}, {3, 4}, {4, 5}, {5, 6}, {6, 7}, {7, 8}, {8, 9}, {9, 10}, {10, 11}, {11, 12}, {12, 13}, {13, 14}, {14, 15}, {15, 16}, {16, 17}, {17, 18}, {18, 19}, {19, 20}, {20, 21}, {21, 22}, {22, 23}, {23, 24}, {24, 28}, {25, 26}, {26, 27}, {27, 25}}
// Hole represented by a point lying inside it
var holes = [][2]float64{
{0.47, -0.5},
}
verts, faces := ConstrainedDelaunay(pts, segs, holes)

const (
numExpectedConstrainedVertices = 29
numExpectedConstrainedFaces = numExpectedConstrainedVertices
)

assertEq(t, numExpectedConstrainedVertices, len(verts),
"ConstrainedDelaunay: Did not get expected number of vertices")
assertEq(t, numExpectedConstrainedFaces, len(faces),
"ConstrainedDelaunay: Did not get expected number of faces")

const (
numExpectedConformingVertices = 61
numExpectedConformingFaces = numExpectedConformingVertices
)

verts, faces = ConformingDelaunay(pts, segs, holes)

assertEq(t, numExpectedConformingVertices, len(verts),
"ConformingDelaunay: Did not get expected number of vertices")
assertEq(t, numExpectedConformingFaces, len(faces),
"ConformingDelaunay: Did not get expected number of faces")

}

func assertEq(t *testing.T, expected, actual int, msg string) {
t.Helper()
if expected != actual {
t.Log(msg)
t.Fatalf("Expected %d, but was %d", expected, actual)
}
}
163 changes: 124 additions & 39 deletions wrapper.go
Original file line number Diff line number Diff line change
Expand Up @@ -6,15 +6,29 @@ package triangle
*/
// #cgo LDFLAGS: -lm
import "C"
import "unsafe"
import (
"unsafe"
)

// The allocations here are using C.malloc and C.free over trimalloc and
// trifree, as C.malloc works the same but gives a better error report in cases
// if we go out of memory.

func trimalloc(size C.ulong) unsafe.Pointer {
return C.malloc(size)
}

func trifree(ptr unsafe.Pointer) {
C.free(ptr)
}

type triangulateIO struct {
ct *C.struct_triangulateio
}

func NewTriangulateIO() *triangulateIO {
t := triangulateIO{}
t.ct = (*C.struct_triangulateio)(C.malloc(C.sizeof_struct_triangulateio))
t.ct = (*C.struct_triangulateio)(trimalloc(C.sizeof_struct_triangulateio))
if t.ct == nil {
panic("Unable to allocate memory")
}
Expand Down Expand Up @@ -45,7 +59,24 @@ func NewTriangulateIO() *triangulateIO {
}

func FreeTriangulateIO(t *triangulateIO) {
C.free(unsafe.Pointer(t.ct))
trifree(unsafe.Pointer(t.ct.edgelist))
trifree(unsafe.Pointer(t.ct.edgemarkerlist))

// trifree(unsafe.Pointer(t.ct.holelist))
// hole list is reused in both in and out of ConstrainedDelaunay

trifree(unsafe.Pointer(t.ct.neighborlist))
trifree(unsafe.Pointer(t.ct.normlist))
trifree(unsafe.Pointer(t.ct.pointattributelist))
trifree(unsafe.Pointer(t.ct.pointlist))
trifree(unsafe.Pointer(t.ct.pointmarkerlist))
trifree(unsafe.Pointer(t.ct.regionlist))
trifree(unsafe.Pointer(t.ct.segmentlist))
trifree(unsafe.Pointer(t.ct.segmentmarkerlist))
trifree(unsafe.Pointer(t.ct.trianglearealist))
trifree(unsafe.Pointer(t.ct.triangleattributelist))
trifree(unsafe.Pointer(t.ct.trianglelist))
trifree(unsafe.Pointer(t.ct))
}

func (t *triangulateIO) NumberOfEdges() int {
Expand All @@ -69,87 +100,93 @@ func (t *triangulateIO) NumberOfTriangles() int {
}

func (t *triangulateIO) Normals() [][2]float64 {
return cArrToFlt64Slice2D(unsafe.Pointer(t.ct.normlist), t.NumberOfEdges())
return cArrToFlt64Slice2D(t.ct.normlist, t.NumberOfEdges())
}

func (t *triangulateIO) Edges() [][2]int32 {
return cArrToIntSlice2D(unsafe.Pointer(t.ct.edgelist), t.NumberOfEdges())
return cArrToIntSlice2D(t.ct.edgelist, t.NumberOfEdges())
}

func (t *triangulateIO) Points() [][2]float64 {
return cArrToFlt64Slice2D(unsafe.Pointer(t.ct.pointlist), t.NumberOfPoints())
return cArrToFlt64Slice2D(t.ct.pointlist, t.NumberOfPoints())
}

func (t *triangulateIO) PointMarkers() []int32 {
return cArrToIntSlice(unsafe.Pointer(t.ct.pointmarkerlist), t.NumberOfPoints())
return cArrToIntSlice(t.ct.pointmarkerlist, t.NumberOfPoints())
}

func (t *triangulateIO) Segments() [][2]int32 {
return cArrToIntSlice2D(unsafe.Pointer(t.ct.segmentlist), t.NumberOfSegments())
return cArrToIntSlice2D(t.ct.segmentlist, t.NumberOfSegments())
}

func (t *triangulateIO) SetEdges(edges [][2]int32) {
t.ct.edgelist = (*C.int)(unsafe.Pointer(&edges[0][0]))
t.ct.numberofedges = C.int(len(edges))
trifree(unsafe.Pointer(t.ct.edgelist))
t.ct.edgelist, t.ct.numberofedges = intSlice2DToCArr(edges)
}

func (t *triangulateIO) SetPoints(pts [][2]float64) {
t.ct.pointlist = (*C.double)(unsafe.Pointer(&pts[0]))
t.ct.numberofpoints = C.int(len(pts))
trifree(unsafe.Pointer(t.ct.pointlist))
t.ct.pointlist, t.ct.numberofpoints = flt64Slice2DToCArr(pts)
}

func (t *triangulateIO) SetPointMarkers(markers []int32) {
t.ct.pointmarkerlist = (*C.int)(unsafe.Pointer(&markers[0]))
trifree(unsafe.Pointer(t.ct.pointmarkerlist))
t.ct.pointmarkerlist, _ = intSliceToCArr(markers)
}

func (t *triangulateIO) SetSegments(segments [][2]int32) {
t.ct.segmentlist = (*C.int)(unsafe.Pointer(&segments[0][0]))
t.ct.numberofsegments = C.int(len(segments))
trifree(unsafe.Pointer(t.ct.segmentlist))
t.ct.segmentlist, t.ct.numberofsegments = intSlice2DToCArr(segments)
}

func (t *triangulateIO) SetSegmentMarkers(markers []int32) {
t.ct.segmentmarkerlist = (*C.int)(unsafe.Pointer(&markers[0]))
trifree(unsafe.Pointer(t.ct.segmentmarkerlist))
t.ct.segmentmarkerlist, _ = intSliceToCArr(markers)
}

func (t *triangulateIO) SetTriangles(tri [][3]int32) {
t.ct.trianglelist = (*C.int)(unsafe.Pointer(&tri[0][0]))
trifree(unsafe.Pointer(t.ct.trianglelist))
t.ct.trianglelist, t.ct.numberoftriangles = intSlice3DToCArr(tri)
}

func (t *triangulateIO) SetTriangleAreas(areas []float64) {
t.ct.trianglearealist = (*C.double)(unsafe.Pointer(&areas[0]))
trifree(unsafe.Pointer(t.ct.trianglearealist))
t.ct.trianglearealist, _ = flt64SliceToCArr(areas)
}

func (t *triangulateIO) SetHoles(holes [][2]float64) {
t.ct.holelist = (*C.double)(unsafe.Pointer(&holes[0][0]))
t.ct.numberofholes = C.int(len(holes))
if t.ct.holelist != nil {
trifree(unsafe.Pointer(t.ct.holelist))
}
t.ct.holelist, t.ct.numberofholes = flt64Slice2DToCArr(holes)
}

func (t *triangulateIO) Triangles() [][3]int32 {
return cArrToIntSlice3D(unsafe.Pointer(t.ct.trianglelist), t.NumberOfTriangles())
return cArrToIntSlice3D(t.ct.trianglelist, t.NumberOfTriangles())
}

func triang(opt string, in, out, vorout *triangulateIO) {
copt := C.CString(opt)
defer C.free(unsafe.Pointer(copt))
defer trifree(unsafe.Pointer(copt))
if vorout == nil {
C.triangulate(copt, in.ct, out.ct, nil)
} else {
C.triangulate(copt, in.ct, out.ct, vorout.ct)
}
}

func cArrToIntSlice(ptr unsafe.Pointer, length int) []int32 {
slice := (*[1 << 30]C.int)(ptr)[:length:length]
func cArrToIntSlice(ptr *C.int, length int) []int32 {
slice := (*[1 << 30]C.int)(unsafe.Pointer(ptr))[:length:length]
result := make([]int32, length)
for i := 0; i < length; i++ {
result[i] = int32(slice[i])
}
return result
}

func cArrToIntSlice2D(ptr unsafe.Pointer, length int) [][2]int32 {
func cArrToIntSlice2D(ptr *C.int, length int) [][2]int32 {
sz := length * 2
slice := (*[1 << 30]C.int)(ptr)[:sz:sz]
slice := (*[1 << 30]C.int)(unsafe.Pointer(ptr))[:sz:sz]
result := make([][2]int32, length)
for i := 0; i < length; i++ {
j := i * 2
Expand All @@ -158,9 +195,9 @@ func cArrToIntSlice2D(ptr unsafe.Pointer, length int) [][2]int32 {
return result
}

func cArrToIntSlice3D(ptr unsafe.Pointer, length int) [][3]int32 {
func cArrToIntSlice3D(ptr *C.int, length int) [][3]int32 {
sz := length * 3
slice := (*[1 << 30]C.int)(ptr)[:sz:sz]
slice := (*[1 << 30]C.int)(unsafe.Pointer(ptr))[:sz:sz]
result := make([][3]int32, length)
for i := 0; i < length; i++ {
j := i * 3
Expand All @@ -169,22 +206,70 @@ func cArrToIntSlice3D(ptr unsafe.Pointer, length int) [][3]int32 {
return result
}

func cArrToFlt64Slice(ptr unsafe.Pointer, length int) []float64 {
slice := (*[1 << 30]C.double)(ptr)[:length:length]
result := make([]float64, length)
for i := 0; i < length; i++ {
result[i] = float64(slice[i])
}
return result
}

func cArrToFlt64Slice2D(ptr unsafe.Pointer, length int) [][2]float64 {
func cArrToFlt64Slice2D(ptr *C.double, length int) [][2]float64 {
sz := length * 2
slice := (*[1 << 30]C.double)(ptr)[:sz:sz]
slice := (*[1 << 30]C.double)(unsafe.Pointer(ptr))[:sz:sz]
result := make([][2]float64, length)
for i := 0; i < length; i++ {
j := i * 2
result[i] = [2]float64{float64(slice[j]), float64(slice[j+1])}
}
return result
}

func intSliceToCArr(slice []int32) (*C.int, C.int) {
sz := len(slice)
ptr := (*C.int)(trimalloc(C.sizeof_int * (C.ulong)(sz)))
cArr := (*[1 << 30]C.int)(unsafe.Pointer(ptr))[:sz:sz]
for i := range slice {
cArr[i] = C.int(slice[i])
}
return ptr, C.int(len(slice))
}

func intSlice2DToCArr(slice [][2]int32) (*C.int, C.int) {
sz := 2 * len(slice)
ptr := (*C.int)(trimalloc(C.sizeof_int * (C.ulong)(sz)))
cArr := (*[1 << 30]C.int)(unsafe.Pointer(ptr))[:sz:sz]
for i := range slice {
j := 2 * i
cArr[j] = C.int(slice[i][0])
cArr[j+1] = C.int(slice[i][1])
}
return ptr, C.int(len(slice))
}

func intSlice3DToCArr(slice [][3]int32) (*C.int, C.int) {
sz := 3 * len(slice)
ptr := (*C.int)(trimalloc(C.sizeof_int * (C.ulong)(sz)))
cArr := (*[1 << 30]C.int)(unsafe.Pointer(ptr))[:sz:sz]
for i := range slice {
j := 3 * i
cArr[j] = C.int(slice[i][0])
cArr[j+1] = C.int(slice[i][1])
cArr[j+2] = C.int(slice[i][2])
}
return ptr, C.int(len(slice))
}

func flt64SliceToCArr(slice []float64) (*C.double, C.int) {
sz := len(slice)
ptr := (*C.double)(trimalloc(C.sizeof_double * (C.ulong)(sz)))
cArr := (*[1 << 30]C.double)(unsafe.Pointer(ptr))[:sz:sz]
for i := range slice {
cArr[i] = C.double(slice[i])
}
return ptr, C.int(len(slice))
}

func flt64Slice2DToCArr(slice [][2]float64) (*C.double, C.int) {
sz := 2 * len(slice)
ptr := (*C.double)(trimalloc(C.sizeof_double * (C.ulong)(sz)))
cArr := (*[1 << 30]C.double)(unsafe.Pointer(ptr))[:sz:sz]
for i := range slice {
j := i * 2
cArr[j] = C.double(slice[i][0])
cArr[j+1] = C.double(slice[i][1])
}
return ptr, C.int(len(slice))
}