forked from scipy/scipy
/
qhull.pyx
1308 lines (1063 loc) · 40.8 KB
/
qhull.pyx
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
"""
Wrappers for Qhull triangulation, plus some additional N-D geometry utilities
.. versionadded:: 0.9
"""
#
# Copyright (C) Pauli Virtanen, 2010.
#
# Distributed under the same BSD license as Scipy.
#
import threading
import numpy as np
cimport numpy as np
cimport cython
cimport qhull
__all__ = ['Delaunay', 'tsearch']
#------------------------------------------------------------------------------
# Qhull interface
#------------------------------------------------------------------------------
cdef extern from "stdio.h":
extern void *stdin
extern void *stderr
extern void *stdout
cdef extern from "math.h":
double fabs(double x) nogil
double sqrt(double x) nogil
cdef extern from "qhull/src/qset.h":
ctypedef union setelemT:
void *p
int i
ctypedef struct setT:
int maxsize
setelemT e[1]
cdef extern from "qhull/src/qhull.h":
ctypedef double realT
ctypedef double coordT
ctypedef double pointT
ctypedef int boolT
ctypedef unsigned int flagT
ctypedef struct facetT:
coordT offset
coordT *center
coordT *normal
facetT *next
facetT *previous
unsigned id
setT *vertices
setT *neighbors
flagT simplicial
flagT flipped
flagT upperdelaunay
ctypedef struct vertexT:
vertexT *next
vertexT *previous
unsigned int id, visitid
pointT *point
setT *neighbours
ctypedef struct qhT:
boolT DELAUNAY
boolT SCALElast
boolT KEEPcoplanar
boolT MERGEexact
boolT NOerrexit
boolT PROJECTdelaunay
boolT ATinfinity
int normal_size
char *qhull_command
facetT *facet_list
facetT *facet_tail
int num_facets
unsigned int facet_id
pointT *first_point
pointT *input_points
realT last_low
realT last_high
realT last_newhigh
realT max_outside
realT MINoutside
realT DISTround
extern qhT qh_qh
extern int qh_PRINToff
extern int qh_ALL
void qh_init_A(void *inp, void *out, void *err, int argc, char **argv) nogil
void qh_init_B(realT *points, int numpoints, int dim, boolT ismalloc) nogil
void qh_checkflags(char *, char *) nogil
void qh_initflags(char *) nogil
void qh_option(char *, char*, char* ) nogil
void qh_freeqhull(boolT) nogil
void qh_memfreeshort(int *curlong, int *totlong) nogil
void qh_qhull() nogil
void qh_check_output() nogil
void qh_produce_output() nogil
void qh_triangulate() nogil
void qh_checkpolygon() nogil
void qh_findgood_all() nogil
void qh_appendprint(int format) nogil
realT *qh_readpoints(int* num, int *dim, boolT* ismalloc) nogil
int qh_new_qhull(int dim, int numpoints, realT *points,
boolT ismalloc, char* qhull_cmd, void *outfile,
void *errfile) nogil
int qh_pointid(pointT *point) nogil
# Qhull is not threadsafe: needs locking
_qhull_lock = threading.Lock()
#------------------------------------------------------------------------------
# LAPACK interface
#------------------------------------------------------------------------------
cdef extern from "qhull_blas.h":
void qh_dgetrf(int *m, int *n, double *a, int *lda, int *ipiv,
int *info) nogil
void qh_dgetrs(char *trans, int *n, int *nrhs, double *a, int *lda,
int *ipiv, double *b, int *ldb, int *info) nogil
void qh_dgecon(char *norm, int *n, double *a, int *lda, double *anorm,
double *rcond, double *work, int *iwork, int *info) nogil
#------------------------------------------------------------------------------
# Delaunay triangulation using Qhull
#------------------------------------------------------------------------------
def _construct_delaunay(np.ndarray[np.double_t, ndim=2] points):
"""
Perform Delaunay triangulation of the given set of points.
"""
# Run qhull with the options
#
# - d: perform delaunay triangulation
# - Qbb: scale last coordinate for Delaunay
# - Qz: reduces Delaunay precision errors for cospherical sites
# - Qt: output only simplical facets (can produce degenerate 0-area ones)
#
cdef char *options = "qhull d Qz Qbb Qt"
cdef int curlong, totlong
cdef int dim
cdef int numpoints
cdef int exitcode
points = np.ascontiguousarray(points)
numpoints = points.shape[0]
dim = points.shape[1]
if numpoints <= 0:
raise ValueError("No points to triangulate")
if dim < 2:
raise ValueError("Need at least 2-D data to triangulate")
_qhull_lock.acquire()
try:
qh_qh.NOerrexit = 1
with nogil:
exitcode = qh_new_qhull(dim, numpoints, <realT*>points.data, 0,
options, NULL, stderr)
try:
if exitcode != 0:
raise RuntimeError("Qhull error")
with nogil:
qh_triangulate() # get rid of non-simplical facets
if qh_qh.SCALElast:
paraboloid_scale = qh_qh.last_newhigh / (
qh_qh.last_high - qh_qh.last_low)
paraboloid_shift = - qh_qh.last_low * paraboloid_scale
else:
paraboloid_scale = 1.0
paraboloid_shift = 0.0
vertices, neighbors, equations = \
_qhull_get_facet_array(dim, numpoints)
return (vertices, neighbors, equations,
paraboloid_scale, paraboloid_shift)
finally:
with nogil:
qh_freeqhull(0)
qh_memfreeshort(&curlong, &totlong)
if curlong != 0 or totlong != 0:
raise RuntimeError("qhull: did not free %d bytes (%d pieces)" %
(totlong, curlong))
finally:
_qhull_lock.release()
@cython.boundscheck(False)
@cython.cdivision(True)
def _qhull_get_facet_array(int ndim, int numpoints):
"""
Return array of simplical facets currently in Qhull.
Returns
-------
vertices : array of int, shape (nfacets, ndim+1)
Indices of coordinates of vertices forming the simplical facets
neighbors : array of int, shape (nfacets, ndim)
Indices of neighboring facets. The kth neighbor is opposite
the kth vertex, and the first neighbor is the horizon facet
for the first vertex.
Facets extending to infinity are denoted with index -1.
"""
cdef facetT* facet
cdef facetT* neighbor
cdef vertexT *vertex
cdef int i, j, point, error_non_simplical
cdef np.ndarray[np.npy_int, ndim=2] vertices
cdef np.ndarray[np.npy_int, ndim=2] neighbors
cdef np.ndarray[np.double_t, ndim=2] equations
cdef np.ndarray[np.npy_int, ndim=1] id_map
id_map = np.empty((qh_qh.facet_id,), dtype=np.intc)
id_map.fill(-1)
# Compute facet indices
facet = qh_qh.facet_list
j = 0
while facet and facet.next:
if facet.simplicial and not facet.upperdelaunay:
id_map[facet.id] = j
j += 1
facet = facet.next
# Allocate output
vertices = np.zeros((j, ndim+1), dtype=np.intc)
neighbors = np.zeros((j, ndim+1), dtype=np.intc)
equations = np.zeros((j, ndim+2), dtype=np.double)
# Retrieve facet information
error_non_simplical = 0
with nogil:
facet = qh_qh.facet_list
j = 0
while facet and facet.next:
if not facet.simplicial:
error_non_simplical = 1
break
if facet.upperdelaunay:
facet = facet.next
continue
# Save vertex info
for i in xrange(ndim+1):
vertex = <vertexT*>facet.vertices.e[i].p
point = qh_pointid(vertex.point)
vertices[j, i] = point
# Save neighbor info
for i in xrange(ndim+1):
neighbor = <facetT*>facet.neighbors.e[i].p
neighbors[j,i] = id_map[neighbor.id]
# Save simplex equation info
for i in xrange(ndim+1):
equations[j,i] = facet.normal[i]
equations[j,ndim+1] = facet.offset
j += 1
facet = facet.next
if error_non_simplical:
raise ValueError("non-simplical facet encountered")
return vertices, neighbors, equations
#------------------------------------------------------------------------------
# Barycentric coordinates
#------------------------------------------------------------------------------
@cython.boundscheck(False)
@cython.cdivision(True)
def _get_barycentric_transforms(np.ndarray[np.double_t, ndim=2] points,
np.ndarray[np.npy_int, ndim=2] vertices,
double eps):
"""
Compute barycentric affine coordinate transformations for given
simplices.
Returns
-------
Tinvs : array, shape (nsimplex, ndim+1, ndim)
Barycentric transforms for each simplex.
Tinvs[i,:ndim,:ndim] contains inverse of the matrix ``T``,
and Tinvs[i,ndim,:] contains the vector ``r_n`` (see below).
Notes
-----
Barycentric transform from ``x`` to ``c`` is defined by::
T c = x - r_n
where the ``r_1, ..., r_n`` are the vertices of the simplex.
The matrix ``T`` is defined by the condition::
T e_j = r_j - r_n
where ``e_j`` is the unit axis vector, e.g, ``e_2 = [0,1,0,0,...]``
This implies that ``T_ij = (r_j - r_n)_i``.
For the barycentric transforms, we need to compute the inverse
matrix ``T^-1`` and store the vectors ``r_n`` for each vertex.
These are stacked into the `Tinvs` returned.
"""
cdef np.ndarray[np.double_t, ndim=2] T
cdef np.ndarray[np.double_t, ndim=3] Tinvs
cdef int isimplex
cdef int i, j, n, nrhs, lda, ldb, info
cdef int ipiv[NPY_MAXDIMS+1]
cdef int ndim, nsimplex
cdef double centroid[NPY_MAXDIMS], c[NPY_MAXDIMS+1]
cdef double *transform
cdef double anorm, rcond
cdef double nan, rcond_limit
cdef double work[4*NPY_MAXDIMS]
cdef int iwork[NPY_MAXDIMS]
cdef double x1, x2, x3
cdef double y1, y2, y3
cdef double det
nan = np.nan
ndim = points.shape[1]
nsimplex = vertices.shape[0]
T = np.zeros((ndim, ndim), dtype=np.double)
Tinvs = np.zeros((nsimplex, ndim+1, ndim), dtype=np.double)
# Maximum inverse condition number to allow: we want at least three
# of the digits be significant, to be safe
rcond_limit = 1000*eps
with nogil:
for isimplex in xrange(nsimplex):
for i in xrange(ndim):
Tinvs[isimplex,ndim,i] = points[vertices[isimplex,ndim],i]
for j in xrange(ndim):
T[i,j] = (points[vertices[isimplex,j],i]
- Tinvs[isimplex,ndim,i])
Tinvs[isimplex,i,i] = 1
# compute 1-norm for estimating condition number
anorm = _matrix_norm1(ndim, <double*>T.data)
# LU decomposition
n = ndim
nrhs = ndim
lda = ndim
ldb = ndim
qh_dgetrf(&n, &n, <double*>T.data, &lda, ipiv, &info)
# Check condition number
if info == 0:
qh_dgecon("1", &n, <double*>T.data, &lda, &anorm, &rcond,
work, iwork, &info)
if rcond < rcond_limit:
# The transform seems singular
info = 1
# Compute transform
if info == 0:
qh_dgetrs("N", &n, &nrhs, <double*>T.data, &lda, ipiv,
(<double*>Tinvs.data) + ndim*(ndim+1)*isimplex,
&ldb, &info)
# Deal with degenerate simplices
if info != 0:
for i in range(ndim+1):
for j in range(ndim):
Tinvs[isimplex,i,j] = nan
return Tinvs
@cython.boundscheck(False)
cdef double _matrix_norm1(int n, double *a) nogil:
"""Compute the 1-norm of a square matrix given in in Fortran order"""
cdef double maxsum = 0, colsum
cdef int i, j
for j in range(n):
colsum = 0
for i in range(n):
colsum += fabs(a[0])
a += 1
if maxsum < colsum:
maxsum = colsum
return maxsum
cdef int _barycentric_inside(int ndim, double *transform,
double *x, double *c, double eps) nogil:
"""
Check whether point is inside a simplex, using barycentric
coordinates. `c` will be filled with barycentric coordinates, if
the point happens to be inside.
"""
cdef int i, j
c[ndim] = 1.0
for i in xrange(ndim):
c[i] = 0
for j in xrange(ndim):
c[i] += transform[ndim*i + j] * (x[j] - transform[ndim*ndim + j])
c[ndim] -= c[i]
if not (-eps <= c[i] <= 1 + eps):
return 0
if not (-eps <= c[ndim] <= 1 + eps):
return 0
return 1
cdef void _barycentric_coordinate_single(int ndim, double *transform,
double *x, double *c, int i) nogil:
"""
Compute a single barycentric coordinate.
Before the ndim+1'th coordinate can be computed, the other must have
been computed earlier.
"""
cdef int j
if i == ndim:
c[ndim] = 1.0
for j in xrange(ndim):
c[ndim] -= c[j]
else:
c[i] = 0
for j in xrange(ndim):
c[i] += transform[ndim*i + j] * (x[j] - transform[ndim*ndim + j])
cdef void _barycentric_coordinates(int ndim, double *transform,
double *x, double *c) nogil:
"""
Compute barycentric coordinates.
"""
cdef int i, j
c[ndim] = 1.0
for i in xrange(ndim):
c[i] = 0
for j in xrange(ndim):
c[i] += transform[ndim*i + j] * (x[j] - transform[ndim*ndim + j])
c[ndim] -= c[i]
#------------------------------------------------------------------------------
# N-D geometry
#------------------------------------------------------------------------------
cdef void _lift_point(DelaunayInfo_t *d, double *x, double *z) nogil:
cdef int i
z[d.ndim] = 0
for i in xrange(d.ndim):
z[i] = x[i]
z[d.ndim] += x[i]**2
z[d.ndim] *= d.paraboloid_scale
z[d.ndim] += d.paraboloid_shift
cdef double _distplane(DelaunayInfo_t *d, int isimplex, double *point) nogil:
"""
qh_distplane
"""
cdef double dist
cdef int k
dist = d.equations[isimplex*(d.ndim+2) + d.ndim+1]
for k in xrange(d.ndim+1):
dist += d.equations[isimplex*(d.ndim+2) + k] * point[k]
return dist
#------------------------------------------------------------------------------
# Iterating over ridges connected to a vertex in 2D
#------------------------------------------------------------------------------
cdef void _RidgeIter2D_init(RidgeIter2D_t *it, DelaunayInfo_t *d,
int vertex) nogil:
"""
Start iteration over all triangles connected to the given vertex.
"""
cdef double c[3]
cdef int k, ivertex, start
start = 0
it.info = d
it.vertex = vertex
it.triangle = d.vertex_to_simplex[vertex]
it.start_triangle = it.triangle
it.restart = 0
if it.triangle != -1:
# find some edge connected to this vertex
for k in xrange(3):
ivertex = it.info.vertices[it.triangle*3 + k]
if ivertex != vertex:
it.vertex2 = ivertex
it.index = k
it.start_index = k
break
else:
it.start_index = -1
it.index = -1
cdef void _RidgeIter2D_next(RidgeIter2D_t *it) nogil:
cdef int itri, k, ivertex
#
# Remember: k-th edge and k-th neigbour are opposite vertex k;
# imagine now we are iterating around vertex `O`
#
# .O------,
# ./ |\. |
# ./ | \. |
# \ | \. |
# \ |k \. |
# \ | \.|
# `+------k
#
if it.restart:
if it.start_index == -1:
# we already did that -> we have iterated over everything
it.index = -1
return
# restart to opposite direction
it.triangle = it.start_triangle
for k in xrange(3):
ivertex = it.info.vertices[it.triangle*3 + k]
if ivertex != it.vertex and k != it.start_index:
it.index = k
it.vertex2 = ivertex
break
it.start_index = -1
it.restart = 0
if it.info.neighbors[it.triangle*3 + it.index] == -1:
it.index = -1
return
else:
_RidgeIter2D_next(it)
if it.index == -1:
return
# jump to the next triangle
itri = it.info.neighbors[it.triangle*3 + it.index]
# if it's outside triangulation, take the last edge, and signal
# restart to the opposite direction
if itri == -1:
for k in xrange(3):
ivertex = it.info.vertices[it.triangle*3 + k]
if ivertex != it.vertex and k != it.index:
it.index = k
it.vertex2 = ivertex
break
it.restart = 1
return
# Find at which index we are now:
#
# it.vertex
# O-------k------.
# | \- /
# | \- E B /
# | \- /
# | A \- /
# +---------´
#
# A = it.triangle
# B = itri
# E = it.index
# O = it.vertex
#
for k in xrange(3):
ivertex = it.info.vertices[itri*3 + k]
if it.info.neighbors[itri*3 + k] != it.triangle and \
ivertex != it.vertex:
it.index = k
it.vertex2 = ivertex
break
it.triangle = itri
# check termination
if it.triangle == it.start_triangle:
it.index = -1
return
cdef class RidgeIter2D(object):
cdef RidgeIter2D_t it
cdef object delaunay
cdef DelaunayInfo_t info
def __init__(self, delaunay, ivertex):
if delaunay.ndim != 2:
raise ValueError("RidgeIter2D supports only 2-D")
self.delaunay = delaunay
_get_delaunay_info(&self.info, delaunay, 0, 1)
_RidgeIter2D_init(&self.it, &self.info, ivertex)
def __iter__(self):
return self
def __next__(self):
if self.it.index == -1:
raise StopIteration()
ret = (self.it.vertex, self.it.vertex2, self.it.index, self.it.triangle)
_RidgeIter2D_next(&self.it)
return ret
#------------------------------------------------------------------------------
# Finding simplices
#------------------------------------------------------------------------------
cdef int _is_point_fully_outside(DelaunayInfo_t *d, double *x,
double eps) nogil:
"""
Is the point outside the bounding box of the triangulation?
"""
cdef int i
for i in xrange(d.ndim):
if x[i] < d.min_bound[i] - eps or x[i] > d.max_bound[i] + eps:
return 1
return 0
cdef int _find_simplex_bruteforce(DelaunayInfo_t *d, double *c,
double *x, double eps,
double eps_broad) nogil:
"""
Find simplex containing point `x` by going through all simplices.
"""
cdef int inside, isimplex
cdef int k, m, ineighbor, iself
cdef double *transform
if _is_point_fully_outside(d, x, eps):
return -1
for isimplex in xrange(d.nsimplex):
transform = d.transform + isimplex*d.ndim*(d.ndim+1)
if transform[0] == transform[0]:
# transform is valid (non-nan)
inside = _barycentric_inside(d.ndim, transform, x, c, eps)
if inside:
return isimplex
else:
# transform is invalid (nan, implying degenerate simplex)
# we replace this inside-check by a check of the neighbors
# with a larger epsilon
for k in xrange(d.ndim+1):
ineighbor = d.neighbors[(d.ndim+1)*isimplex + k]
if ineighbor == -1:
continue
transform = d.transform + ineighbor*d.ndim*(d.ndim+1)
if transform[0] != transform[0]:
# another bad simplex
continue
_barycentric_coordinates(d.ndim, transform, x, c)
# Check that the point lies (almost) inside the
# neigbor simplex
inside = 1
for m in xrange(d.ndim+1):
if d.neighbors[(d.ndim+1)*ineighbor + m] == isimplex:
# allow extra leeway towards isimplex
if not (-eps_broad <= c[m] <= 1 + eps):
inside = 0
break
else:
# normal check
if not (-eps <= c[m] <= 1 + eps):
inside = 0
break
if inside:
return ineighbor
return -1
cdef int _find_simplex_directed(DelaunayInfo_t *d, double *c,
double *x, int *start, double eps,
double eps_broad) nogil:
"""
Find simplex containing point `x` via a directed walk in the tesselation.
If the simplex is found, the array `c` is filled with the corresponding
barycentric coordinates.
Notes
-----
The idea here is the following:
1) In a simplex, the k-th neighbour is opposite the k-th vertex.
Call the ridge between them the k-th ridge.
2) If the k-th barycentric coordinate of the target point is negative,
then the k-th vertex and the target point lie on the opposite sides
of the k-th ridge.
3) Consequently, the k-th neighbour simplex is *closer* to the target point
than the present simplex, if projected on the normal of the k-th ridge.
4) In a regular tesselation, hopping to any such direction is OK.
Also, if one of the negative-coordinate neighbors happens to be -1,
then the target point is outside the tesselation (because the
tesselation is convex!).
5) If all barycentric coordinates are in [-eps, 1+eps], we have found the
simplex containing the target point.
6) If all barycentric coordinates are non-negative but 5) is not true,
we are in an inconsistent situation -- this should never happen.
This may however enter an infinite loop due to rounding errors in
the computation of the barycentric coordinates, so the iteration
count needs to be limited, and a fallback to brute force provided.
"""
cdef int k, m, ndim, inside, isimplex, cycle_k
cdef double *transform
cdef double v
ndim = d.ndim
isimplex = start[0]
if isimplex < 0 or isimplex >= d.nsimplex:
isimplex = 0
# The maximum iteration count: it should be large enough so that
# the algorithm usually succeeds, but smaller than nsimplex so
# that for the cases where the algorithm fails, the main cost
# still comes from the brute force search.
for cycle_k in range(1 + d.nsimplex//4):
if isimplex == -1:
break
transform = d.transform + isimplex*ndim*(ndim+1)
inside = 1
for k in xrange(ndim+1):
_barycentric_coordinate_single(ndim, transform, x, c, k)
if c[k] < -eps:
# The target point is in the direction of neighbor `k`!
m = d.neighbors[(ndim+1)*isimplex + k]
if m == -1:
# The point is outside the triangulation: bail out
start[0] = isimplex
return -1
isimplex = m
inside = -1
break
elif c[k] <= 1 + eps:
# we're inside this simplex
pass
else:
# we're outside (or the coordinate is nan; a degenerate simplex)
inside = 0
if inside == -1:
# hopped to another simplex
continue
elif inside == 1:
# we've found the right one!
break
else:
# we've failed utterly (degenerate simplices in the way).
# fall back to brute force
isimplex = _find_simplex_bruteforce(d, c, x, eps, eps_broad)
break
else:
# the algorithm failed to converge -- fall back to brute force
isimplex = _find_simplex_bruteforce(d, c, x, eps, eps_broad)
start[0] = isimplex
return isimplex
cdef int _find_simplex(DelaunayInfo_t *d, double *c,
double *x, int *start, double eps,
double eps_broad) nogil:
"""
Find simplex containing point `x` by walking the triangulation.
Notes
-----
This algorithm is similar as used by ``qh_findbest``. The idea
is the following:
1. Delaunay triangulation is a projection of the lower half of a convex
hull, of points lifted on a paraboloid.
Simplices in the triangulation == facets on the convex hull.
2. If a point belongs to a given simplex in the triangulation,
its image on the paraboloid is on the positive side of
the corresponding facet.
3. However, it is not necessarily the *only* such facet.
4. Also, it is not necessarily the facet whose hyperplane distance
to the point on the paraboloid is the largest.
..note::
If I'm not mistaken, `qh_findbestfacet` finds a facet for
which the plane distance is maximized -- so it doesn't always
return the simplex containing the point given. For example:
>>> p = np.array([(1 - 1e-4, 0.1)])
>>> points = np.array([(0,0), (1, 1), (1, 0), (0.99189033, 0.37674127),
... (0.99440079, 0.45182168)], dtype=np.double)
>>> tri = qhull.delaunay(points)
>>> tri.vertices
array([[4, 1, 0],
[4, 2, 1],
[3, 2, 0],
[3, 4, 0],
[3, 4, 2]])
>>> dist = qhull.plane_distance(tri, p)
>>> dist
array([[-0.12231439, 0.00184863, 0.01049659, -0.04714842,
0.00425905]])
>>> tri.vertices[dist.argmax()]
array([3, 2, 0]
Now, the maximally positive-distant simplex is [3, 2, 0], although
the simplex containing the point is [4, 2, 1].
In this algorithm, we walk around the tesselation trying to locate
a positive-distant facet. After finding one, we fall back to a
directed search.
"""
cdef int isimplex, i, j, k, inside, ineigh, neighbor_found
cdef int ndim
cdef double z[NPY_MAXDIMS+1]
cdef double best_dist, dist
cdef int changed
if _is_point_fully_outside(d, x, eps):
return -1
if d.nsimplex <= 0:
return -1
ndim = d.ndim
isimplex = start[0]
if isimplex < 0 or isimplex >= d.nsimplex:
isimplex = 0
# Lift point to paraboloid
_lift_point(d, x, z)
# Walk the tesselation searching for a facet with a positive planar distance
best_dist = _distplane(d, isimplex, z)
changed = 1
while changed:
if best_dist > 0:
break
changed = 0
for k in xrange(ndim+1):
ineigh = d.neighbors[(ndim+1)*isimplex + k]
if ineigh == -1:
continue
dist = _distplane(d, ineigh, z)
# Note addition of eps -- otherwise, this code does not
# necessarily terminate! The compiler may use extended
# accuracy of the FPU so that (dist > best_dist), but
# after storing to double size, dist == best_dist,
# resulting to non-terminating loop
if dist > best_dist + eps*(1 + fabs(best_dist)):
# Note: this is intentional: we jump in the middle of the cycle,
# and continue the cycle from the next k.
#
# This apparently sweeps the different directions more
# efficiently. We don't need full accuracy, since we do
# a directed search afterwards in any case.
isimplex = ineigh
best_dist = dist
changed = 1
# We should now be somewhere near the simplex containing the point,
# locate it with a directed search
start[0] = isimplex
return _find_simplex_directed(d, c, x, start, eps, eps_broad)
#------------------------------------------------------------------------------
# Delaunay triangulation interface, for Python
#------------------------------------------------------------------------------
class Delaunay(object):
"""
Delaunay(points)
Delaunay tesselation in N dimensions.
Parameters
----------
points : ndarray of floats, shape (npoints, ndim)
Coordinates of points to triangulate
Attributes
----------
points : ndarray of double, shape (npoints, ndim)
Points in the triangulation.
vertices : ndarray of ints, shape (nsimplex, ndim+1)
Indices of vertices forming simplices in the triangulation.
neighbors : ndarray of ints, shape (nsimplex, ndim+1)
Indices of neighbor simplices for each simplex.
The kth neighbor is opposite to the kth vertex.
For simplices at the boundary, -1 denotes no neighbor.
equations : ndarray of double, shape (nsimplex, ndim+2)
[normal, offset] forming the hyperplane equation of the facet
on the paraboloid (see [Qhull]_ documentation for more).
paraboloid_scale, paraboloid_shift : float
Scale and shift for the extra paraboloid dimension
(see [Qhull]_ documentation for more).
transform : ndarray of double, shape (nsimplex, ndim+1, ndim)
Affine transform from ``x`` to the barycentric coordinates ``c``.
This is defined by::
T c = x - r
At vertex ``j``, ``c_j = 1`` and the other coordinates zero.
For simplex ``i``, ``transform[i,:ndim,:ndim]`` contains
inverse of the matrix ``T``, and ``transform[i,ndim,:]``
contains the vector ``r``.
vertex_to_simplex : ndarray of int, shape (npoints,)
Lookup array, from a vertex, to some simplex which it is a part of.
convex_hull : ndarray of int, shape (nfaces, ndim)
Vertices of facets forming the convex hull of the point set.
The array contains the indices of the points belonging to
the (N-1)-dimensional facets that form the convex hull
of the triangulation.
Notes
-----
The tesselation is computed using the Qhull libary [Qhull]_.
.. versionadded:: 0.9
References
----------
.. [Qhull] http://www.qhull.org/
"""
def __init__(self, points):
points = np.ascontiguousarray(points).astype(np.double)
vertices, neighbors, equations, paraboloid_scale, paraboloid_shift = \
_construct_delaunay(points)
self.ndim = points.shape[1]
self.npoints = points.shape[0]
self.nsimplex = vertices.shape[0]