-
Notifications
You must be signed in to change notification settings - Fork 2
/
stripack.f90
7496 lines (6909 loc) · 198 KB
/
stripack.f90
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
subroutine addnod ( nst, k, x, y, z, list, lptr, lend, lnew, ier )
!*****************************************************************************80
!
!! ADDNOD adds a node to a triangulation.
!
! Discussion:
!
! This subroutine adds node K to a triangulation of the
! convex hull of nodes 1, ..., K-1, producing a triangulation
! of the convex hull of nodes 1, ..., K.
!
! The algorithm consists of the following steps: node K
! is located relative to the triangulation (TRFIND), its
! index is added to the data structure (INTADD or BDYADD),
! and a sequence of swaps (SWPTST and SWAP) are applied to
! the arcs opposite K so that all arcs incident on node K
! and opposite node K are locally optimal (satisfy the circumcircle test).
!
! Thus, if a Delaunay triangulation of nodes 1 through K-1 is input,
! a Delaunay triangulation of nodes 1 through K will be output.
!
! Modified:
!
! 15 May 2007
!
! Author:
!
! Robert Renka
!
! Reference:
!
! Robert Renka,
! Algorithm 772: STRIPACK,
! Delaunay Triangulation and Voronoi Diagram on the Surface of a Sphere,
! ACM Transactions on Mathematical Software,
! Volume 23, Number 3, September 1997, pages 416-434.
!
! Parameters:
!
! Input, integer ( kind = 4 ) NST, the index of a node at which TRFIND
! begins its search. Search time depends on the proximity of this node to
! K. If NST < 1, the search is begun at node K-1.
!
! Input, integer ( kind = 4 ) K, the nodal index (index for X, Y, Z, and
! LEND) of the new node to be added. 4 <= K.
!
! Input, real ( kind = 8 ) X(K), Y(K), Z(K), the coordinates of the nodes.
!
! Input/output, integer ( kind = 4 ) LIST(6*(N-2)), LPTR(6*(N-2)), LEND(K),
! LNEW. On input, the data structure associated with the triangulation of
! nodes 1 to K-1. On output, the data has been updated to include node
! K. The array lengths are assumed to be large enough to add node K.
! Refer to TRMESH.
!
! Output, integer ( kind = 4 ) IER, error indicator:
! 0 if no errors were encountered.
! -1 if K is outside its valid range on input.
! -2 if all nodes (including K) are collinear (lie on a common geodesic).
! L if nodes L and K coincide for some L < K.
!
! Local parameters:
!
! B1,B2,B3 = Unnormalized barycentric coordinates returned by TRFIND.
! I1,I2,I3 = Vertex indexes of a triangle containing K
! IN1 = Vertex opposite K: first neighbor of IO2
! that precedes IO1. IN1,IO1,IO2 are in
! counterclockwise order.
! IO1,IO2 = Adjacent neighbors of K defining an arc to
! be tested for a swap
! IST = Index of node at which TRFIND begins its search
! KK = Local copy of K
! KM1 = K-1
! L = Vertex index (I1, I2, or I3) returned in IER
! if node K coincides with a vertex
! LP = LIST pointer
! LPF = LIST pointer to the first neighbor of K
! LPO1 = LIST pointer to IO1
! LPO1S = Saved value of LPO1
! P = Cartesian coordinates of node K
!
implicit none
integer ( kind = 4 ) k
real ( kind = 8 ) b1
real ( kind = 8 ) b2
real ( kind = 8 ) b3
integer ( kind = 4 ) i1
integer ( kind = 4 ) i2
integer ( kind = 4 ) i3
integer ( kind = 4 ) ier
integer ( kind = 4 ) in1
integer ( kind = 4 ) io1
integer ( kind = 4 ) io2
integer ( kind = 4 ) ist
integer ( kind = 4 ) kk
integer ( kind = 4 ) km1
integer ( kind = 4 ) l
integer ( kind = 4 ) lend(k)
integer ( kind = 4 ) list(*)
integer ( kind = 4 ) lnew
integer ( kind = 4 ) lp
integer ( kind = 4 ) lpf
integer ( kind = 4 ) lpo1
integer ( kind = 4 ) lpo1s
integer ( kind = 4 ) lptr(*)
integer ( kind = 4 ) lstptr
integer ( kind = 4 ) nst
real ( kind = 8 ) p(3)
logical swptst
real ( kind = 8 ) x(k)
real ( kind = 8 ) y(k)
real ( kind = 8 ) z(k)
kk = k
if ( kk < 4 ) then
ier = -1
write ( *, '(a)' ) ' '
write ( *, '(a)' ) 'ADDNOD - Fatal error!'
write ( *, '(a)' ) ' K < 4.'
stop
end if
!
! Initialization:
!
km1 = kk - 1
ist = nst
if ( ist < 1 ) then
ist = km1
end if
p(1) = x(kk)
p(2) = y(kk)
p(3) = z(kk)
!
! Find a triangle (I1,I2,I3) containing K or the rightmost
! (I1) and leftmost (I2) visible boundary nodes as viewed
! from node K.
!
call trfind ( ist, p, km1, x, y, z, list, lptr, lend, b1, b2, b3, &
i1, i2, i3 )
!
! Test for collinear or duplicate nodes.
!
if ( i1 == 0 ) then
ier = -2
write ( *, '(a)' ) ' '
write ( *, '(a)' ) 'ADDNOD - Fatal error!'
write ( *, '(a)' ) ' The nodes are coplanar.'
stop
end if
if ( i3 /= 0 ) then
l = i1
if ( p(1) == x(l) .and. p(2) == y(l) .and. p(3) == z(l) ) then
ier = l
write ( *, '(a)' ) ' '
write ( *, '(a)' ) 'ADDNOD - Fatal error!'
write ( *, '(a,i8,a,i8)' ) ' Node ', l, ' is equal to node ', k
stop
end if
l = i2
if ( p(1) == x(l) .and. p(2) == y(l) .and. p(3) == z(l) ) then
ier = l
write ( *, '(a)' ) ' '
write ( *, '(a)' ) 'ADDNOD - Fatal error!'
write ( *, '(a,i8,a,i8)' ) ' Node ', l, ' is equal to node ', k
stop
end if
l = i3
if ( p(1) == x(l) .and. p(2) == y(l) .and. p(3) == z(l) ) then
ier = l
write ( *, '(a)' ) ' '
write ( *, '(a)' ) 'ADDNOD - Fatal error!'
write ( *, '(a,i8,a,i8)' ) ' Node ', l, ' is equal to node ', k
stop
end if
call intadd ( kk, i1, i2, i3, list, lptr, lend, lnew )
else
if ( i1 /= i2 ) then
call bdyadd ( kk, i1,i2, list, lptr, lend, lnew )
else
call covsph ( kk, i1, list, lptr, lend, lnew )
end if
end if
ier = 0
!
! Initialize variables for optimization of the triangulation.
!
lp = lend(kk)
lpf = lptr(lp)
io2 = list(lpf)
lpo1 = lptr(lpf)
io1 = abs ( list(lpo1) )
!
! Begin loop: find the node opposite K.
!
do
lp = lstptr ( lend(io1), io2, list, lptr )
if ( 0 <= list(lp) ) then
lp = lptr(lp)
in1 = abs ( list(lp) )
!
! Swap test: if a swap occurs, two new arcs are
! opposite K and must be tested.
!
lpo1s = lpo1
if ( .not. swptst ( in1, kk, io1, io2, x, y, z ) ) then
if ( lpo1 == lpf .or. list(lpo1) < 0 ) then
exit
end if
io2 = io1
lpo1 = lptr(lpo1)
io1 = abs ( list(lpo1) )
cycle
end if
call swap ( in1, kk, io1, io2, list, lptr, lend, lpo1 )
!
! A swap is not possible because KK and IN1 are already
! adjacent. This error in SWPTST only occurs in the
! neutral case and when there are nearly duplicate nodes.
!
if ( lpo1 /= 0 ) then
io1 = in1
cycle
end if
lpo1 = lpo1s
end if
!
! No swap occurred. Test for termination and reset IO2 and IO1.
!
if ( lpo1 == lpf .or. list(lpo1) < 0 ) then
exit
end if
io2 = io1
lpo1 = lptr(lpo1)
io1 = abs ( list(lpo1) )
end do
return
end
function arc_cosine ( c )
!*****************************************************************************80
!
!! ARC_COSINE computes the arc cosine function, with argument truncation.
!
! Discussion:
!
! If you call your system ACOS routine with an input argument that is
! outside the range [-1.0, 1.0 ], you may get an unpleasant surprise.
! This routine truncates arguments outside the range.
!
! Licensing:
!
! This code is distributed under the GNU LGPL license.
!
! Modified:
!
! 02 December 2000
!
! Author:
!
! John Burkardt
!
! Parameters:
!
! Input, real ( kind = 8 ) C, the argument.
!
! Output, real ( kind = 8 ) ARC_COSINE, an angle whose cosine is C.
!
implicit none
real ( kind = 8 ) arc_cosine
real ( kind = 8 ) c
real ( kind = 8 ) c2
c2 = c
c2 = max ( c2, -1.0D+00 )
c2 = min ( c2, +1.0D+00 )
arc_cosine = acos ( c2 )
return
end
function areas ( v1, v2, v3 )
!*****************************************************************************80
!
!! AREAS computes the area of a spherical triangle on the unit sphere.
!
! Discussion:
!
! This function returns the area of a spherical triangle
! on the unit sphere.
!
! Modified:
!
! 16 June 2007
!
! Author:
!
! Robert Renka
!
! Reference:
!
! Robert Renka,
! Algorithm 772: STRIPACK,
! Delaunay Triangulation and Voronoi Diagram on the Surface of a Sphere,
! ACM Transactions on Mathematical Software,
! Volume 23, Number 3, September 1997, pages 416-434.
!
! Parameters:
!
! Input, real ( kind = 8 ) V1(3), V2(3), V3(3), the Cartesian coordinates
! of unit vectors (the three triangle vertices in any order). These
! vectors, if nonzero, are implicitly scaled to have length 1.
!
! Output, real ( kind = 8 ) AREAS, the area of the spherical triangle
! defined by V1, V2, and V3, in the range 0 to 2*PI (the area of a
! hemisphere). AREAS = 0 (or 2*PI) if and only if V1, V2, and V3 lie in (or
! close to) a plane containing the origin.
!
! Local parameters:
!
! A1,A2,A3 = Interior angles of the spherical triangle.
!
! CA1,CA2,CA3 = cos(A1), cos(A2), and cos(A3), respectively.
!
! DV1,DV2,DV3 = copies of V1, V2, and V3.
!
! I = DO-loop index and index for Uij.
!
! S12,S23,S31 = Sum of squared components of U12, U23, U31.
!
! U12,U23,U31 = Unit normal vectors to the planes defined by
! pairs of triangle vertices.
!
implicit none
real ( kind = 8 ) a1
real ( kind = 8 ) a2
real ( kind = 8 ) a3
real ( kind = 8 ) areas
real ( kind = 8 ) ca1
real ( kind = 8 ) ca2
real ( kind = 8 ) ca3
real ( kind = 8 ) dv1(3)
real ( kind = 8 ) dv2(3)
real ( kind = 8 ) dv3(3)
real ( kind = 8 ) s12
real ( kind = 8 ) s23
real ( kind = 8 ) s31
real ( kind = 8 ) u12(3)
real ( kind = 8 ) u23(3)
real ( kind = 8 ) u31(3)
real ( kind = 8 ) v1(3)
real ( kind = 8 ) v2(3)
real ( kind = 8 ) v3(3)
dv1(1:3) = v1(1:3)
dv2(1:3) = v2(1:3)
dv3(1:3) = v3(1:3)
!
! Compute cross products Uij = Vi X Vj.
!
u12(1) = dv1(2) * dv2(3) - dv1(3) * dv2(2)
u12(2) = dv1(3) * dv2(1) - dv1(1) * dv2(3)
u12(3) = dv1(1) * dv2(2) - dv1(2) * dv2(1)
u23(1) = dv2(2) * dv3(3) - dv2(3) * dv3(2)
u23(2) = dv2(3) * dv3(1) - dv2(1) * dv3(3)
u23(3) = dv2(1) * dv3(2) - dv2(2) * dv3(1)
u31(1) = dv3(2) * dv1(3) - dv3(3) * dv1(2)
u31(2) = dv3(3) * dv1(1) - dv3(1) * dv1(3)
u31(3) = dv3(1) * dv1(2) - dv3(2) * dv1(1)
!
! Normalize Uij to unit vectors.
!
s12 = dot_product ( u12(1:3), u12(1:3) )
s23 = dot_product ( u23(1:3), u23(1:3) )
s31 = dot_product ( u31(1:3), u31(1:3) )
!
! Test for a degenerate triangle associated with collinear vertices.
!
if ( s12 == 0.0D+00 .or. s23 == 0.0D+00 .or. s31 == 0.0D+00 ) then
areas = 0.0D+00
return
end if
s12 = sqrt ( s12 )
s23 = sqrt ( s23 )
s31 = sqrt ( s31 )
u12(1:3) = u12(1:3) / s12
u23(1:3) = u23(1:3) / s23
u31(1:3) = u31(1:3) / s31
!
! Compute interior angles Ai as the dihedral angles between planes:
! CA1 = cos(A1) = -<U12,U31>
! CA2 = cos(A2) = -<U23,U12>
! CA3 = cos(A3) = -<U31,U23>
!
ca1 = - dot_product ( u12(1:3), u31(1:3) )
ca2 = - dot_product ( u23(1:3), u12(1:3) )
ca3 = - dot_product ( u31(1:3), u23(1:3) )
ca1 = max ( ca1, -1.0D+00 )
ca1 = min ( ca1, +1.0D+00 )
ca2 = max ( ca2, -1.0D+00 )
ca2 = min ( ca2, +1.0D+00 )
ca3 = max ( ca3, -1.0D+00 )
ca3 = min ( ca3, +1.0D+00 )
a1 = acos ( ca1 )
a2 = acos ( ca2 )
a3 = acos ( ca3 )
!
! Compute AREAS = A1 + A2 + A3 - PI.
!
areas = a1 + a2 + a3 - acos ( -1.0D+00 )
if ( areas < 0.0D+00 ) then
areas = 0.0D+00
end if
return
end
subroutine bdyadd ( kk, i1, i2, list, lptr, lend, lnew )
!*****************************************************************************80
!
!! BDYADD adds a boundary node to a triangulation.
!
! Discussion:
!
! This subroutine adds a boundary node to a triangulation
! of a set of KK-1 points on the unit sphere. The data
! structure is updated with the insertion of node KK, but no
! optimization is performed.
!
! This routine is identical to the similarly named routine
! in TRIPACK.
!
! Modified:
!
! 16 June 2007
!
! Author:
!
! Robert Renka
!
! Reference:
!
! Robert Renka,
! Algorithm 772: STRIPACK,
! Delaunay Triangulation and Voronoi Diagram on the Surface of a Sphere,
! ACM Transactions on Mathematical Software,
! Volume 23, Number 3, September 1997, pages 416-434.
!
! Parameters:
!
! Input, integer ( kind = 4 ) KK, the index of a node to be connected to
! the sequence of all visible boundary nodes. 1 <= KK and
! KK must not be equal to I1 or I2.
!
! Input, integer ( kind = 4 ) I1, the first (rightmost as viewed from KK)
! boundary node in the triangulation that is visible from
! node KK (the line segment KK-I1 intersects no arcs.
!
! Input, integer ( kind = 4 ) I2, the last (leftmost) boundary node that
! is visible from node KK. I1 and I2 may be determined by TRFIND.
!
! Input/output, integer ( kind = 4 ) LIST(6*(N-2)), LPTR(6*(N-2)), LEND(N),
! LNEW, the triangulation data structure created by TRMESH.
! Nodes I1 and I2 must be included
! in the triangulation. On output, the data structure is updated with
! the addition of node KK. Node KK is connected to I1, I2, and
! all boundary nodes in between.
!
! Local parameters:
!
! K = Local copy of KK
! LP = LIST pointer
! LSAV = LIST pointer
! N1,N2 = Local copies of I1 and I2, respectively
! NEXT = Boundary node visible from K
! NSAV = Boundary node visible from K
!
implicit none
integer ( kind = 4 ) i1
integer ( kind = 4 ) i2
integer ( kind = 4 ) k
integer ( kind = 4 ) kk
integer ( kind = 4 ) lend(*)
integer ( kind = 4 ) list(*)
integer ( kind = 4 ) lnew
integer ( kind = 4 ) lp
integer ( kind = 4 ) lptr(*)
integer ( kind = 4 ) lsav
integer ( kind = 4 ) n1
integer ( kind = 4 ) n2
integer ( kind = 4 ) next
integer ( kind = 4 ) nsav
k = kk
n1 = i1
n2 = i2
!
! Add K as the last neighbor of N1.
!
lp = lend(n1)
lsav = lptr(lp)
lptr(lp) = lnew
list(lnew) = -k
lptr(lnew) = lsav
lend(n1) = lnew
lnew = lnew + 1
next = -list(lp)
list(lp) = next
nsav = next
!
! Loop on the remaining boundary nodes between N1 and N2,
! adding K as the first neighbor.
!
do
lp = lend(next)
call insert ( k, lp, list, lptr, lnew )
if ( next == n2 ) then
exit
end if
next = -list(lp)
list(lp) = next
end do
!
! Add the boundary nodes between N1 and N2 as neighbors of node K.
!
lsav = lnew
list(lnew) = n1
lptr(lnew) = lnew + 1
lnew = lnew + 1
next = nsav
do
if ( next == n2 ) then
exit
end if
list(lnew) = next
lptr(lnew) = lnew + 1
lnew = lnew + 1
lp = lend(next)
next = list(lp)
end do
list(lnew) = -n2
lptr(lnew) = lsav
lend(k) = lnew
lnew = lnew + 1
return
end
subroutine bnodes ( n, list, lptr, lend, nodes, nb, na, nt )
!*****************************************************************************80
!
!! BNODES returns the boundary nodes of a triangulation.
!
! Discussion:
!
! Given a triangulation of N nodes on the unit sphere created by TRMESH,
! this subroutine returns an array containing the indexes (if any) of
! the counterclockwise sequence of boundary nodes, that is, the nodes on
! the boundary of the convex hull of the set of nodes. The
! boundary is empty if the nodes do not lie in a single
! hemisphere. The numbers of boundary nodes, arcs, and
! triangles are also returned.
!
! Modified:
!
! 16 June 2007
!
! Author:
!
! Robert Renka
!
! Reference:
!
! Robert Renka,
! Algorithm 772: STRIPACK,
! Delaunay Triangulation and Voronoi Diagram on the Surface of a Sphere,
! ACM Transactions on Mathematical Software,
! Volume 23, Number 3, September 1997, pages 416-434.
!
! Parameters:
!
! Input, integer ( kind = 4 ) N, the number of nodes in the triangulation.
! 3 <= N.
!
! Input, integer ( kind = 4 ) LIST(6*(N-2)), LPTR(6*(N-2)), LEND(N), the
! data structure defining the triangulation, created by TRMESH.
!
! Output, integer ( kind = 4 ) NODES(*), the ordered sequence of NB boundary
! node indexes in the range 1 to N. For safety, the dimension of NODES
! should be N.
!
! Output, integer ( kind = 4 ) NB, the number of boundary nodes.
!
! Output, integer ( kind = 4 ) NA, NT, the number of arcs and triangles,
! respectively, in the triangulation.
!
! Local parameters:
!
! K = NODES index
! LP = LIST pointer
! N0 = Boundary node to be added to NODES
! NN = Local copy of N
! NST = First element of nodes (arbitrarily chosen to be
! the one with smallest index)
!
implicit none
integer ( kind = 4 ) n
integer ( kind = 4 ) i
integer ( kind = 4 ) k
integer ( kind = 4 ) lend(n)
integer ( kind = 4 ) list(6*(n-2))
integer ( kind = 4 ) lp
integer ( kind = 4 ) lptr(6*(n-2))
integer ( kind = 4 ) n0
integer ( kind = 4 ) na
integer ( kind = 4 ) nb
integer ( kind = 4 ) nn
integer ( kind = 4 ) nodes(*)
integer ( kind = 4 ) nst
integer ( kind = 4 ) nt
nn = n
!
! Search for a boundary node.
!
nst = 0
do i = 1, nn
lp = lend(i)
if ( list(lp) < 0 ) then
nst = i
exit
end if
end do
!
! The triangulation contains no boundary nodes.
!
if ( nst == 0 ) then
nb = 0
na = 3 * ( nn - 2 )
nt = 2 * ( nn - 2 )
return
end if
!
! NST is the first boundary node encountered.
!
! Initialize for traversal of the boundary.
!
nodes(1) = nst
k = 1
n0 = nst
!
! Traverse the boundary in counterclockwise order.
!
do
lp = lend(n0)
lp = lptr(lp)
n0 = list(lp)
if ( n0 == nst ) then
exit
end if
k = k + 1
nodes(k) = n0
end do
!
! Store the counts.
!
nb = k
nt = 2 * n - nb - 2
na = nt + n - 1
return
end
subroutine circum ( v1, v2, v3, c, ier )
!*****************************************************************************80
!
!! CIRCUM returns the circumcenter of a spherical triangle.
!
! Discussion:
!
! This subroutine returns the circumcenter of a spherical triangle on the
! unit sphere: the point on the sphere surface that is equally distant
! from the three triangle vertices and lies in the same hemisphere, where
! distance is taken to be arc-length on the sphere surface.
!
! Modified:
!
! 16 June 2007
!
! Author:
!
! Robert Renka
!
! Reference:
!
! Robert Renka,
! Algorithm 772: STRIPACK,
! Delaunay Triangulation and Voronoi Diagram on the Surface of a Sphere,
! ACM Transactions on Mathematical Software,
! Volume 23, Number 3, September 1997, pages 416-434.
!
! Parameters:
!
! Input, real ( kind = 8 ) V1(3), V2(3), V3(3), the coordinates of the
! three triangle vertices (unit vectors) in counter clockwise order.
!
! Output, real ( kind = 8 ) C(3), the coordinates of the circumcenter unless
! 0 < IER, in which case C is not defined. C = (V2-V1) X (V3-V1)
! normalized to a unit vector.
!
! Output, integer ( kind = 4 ) IER = Error indicator:
! 0, if no errors were encountered.
! 1, if V1, V2, and V3 lie on a common line: (V2-V1) X (V3-V1) = 0.
!
! Local parameters:
!
! CNORM = Norm of CU: used to compute C
! CU = Scalar multiple of C: E1 X E2
! E1,E2 = Edges of the underlying planar triangle:
! V2-V1 and V3-V1, respectively
! I = DO-loop index
!
implicit none
real ( kind = 8 ) c(3)
real ( kind = 8 ) cnorm
real ( kind = 8 ) cu(3)
real ( kind = 8 ) e1(3)
real ( kind = 8 ) e2(3)
integer ( kind = 4 ) ier
real ( kind = 8 ) v1(3)
real ( kind = 8 ) v2(3)
real ( kind = 8 ) v3(3)
ier = 0
e1(1:3) = v2(1:3) - v1(1:3)
e2(1:3) = v3(1:3) - v1(1:3)
!
! Compute CU = E1 X E2 and CNORM**2.
!
cu(1) = e1(2) * e2(3) - e1(3) * e2(2)
cu(2) = e1(3) * e2(1) - e1(1) * e2(3)
cu(3) = e1(1) * e2(2) - e1(2) * e2(1)
cnorm = sqrt ( sum ( cu(1:3)**2 ) )
!
! The vertices lie on a common line if and only if CU is the zero vector.
!
if ( cnorm == 0.0D+00 ) then
ier = 1
return
end if
c(1:3) = cu(1:3) / cnorm
return
end
subroutine covsph ( kk, n0, list, lptr, lend, lnew )
!*****************************************************************************80
!
!! COVSPH connects an exterior node to boundary nodes, covering the sphere.
!
! Discussion:
!
! This subroutine connects an exterior node KK to all
! boundary nodes of a triangulation of KK-1 points on the
! unit sphere, producing a triangulation that covers the
! sphere. The data structure is updated with the addition
! of node KK, but no optimization is performed. All
! boundary nodes must be visible from node KK.
!
! Modified:
!
! 16 June 2007
!
! Author:
!
! Robert Renka
!
! Reference:
!
! Robert Renka,
! Algorithm 772: STRIPACK,
! Delaunay Triangulation and Voronoi Diagram on the Surface of a Sphere,
! ACM Transactions on Mathematical Software,
! Volume 23, Number 3, September 1997, pages 416-434.
!
! Parameters:
!
! Input, integer ( kind = 4 ) KK = Index of the node to be connected to the
! set of all boundary nodes. 4 <= KK.
!
! Input, integer ( kind = 4 ) N0 = Index of a boundary node (in the range
! 1 to KK-1). N0 may be determined by TRFIND.
!
! Input/output, integer ( kind = 4 ) LIST(6*(N-2)), LPTR(6*(N-2)), LEND(N),
! LNEW, the triangulation data structure created by TRMESH. Node N0 must
! be included in the triangulation. On output, updated with the addition
! of node KK as the last entry. The updated triangulation contains no
! boundary nodes.
!
! Local parameters:
!
! K = Local copy of KK
! LP = LIST pointer
! LSAV = LIST pointer
! NEXT = Boundary node visible from K
! NST = Local copy of N0
!
implicit none
integer ( kind = 4 ) k
integer ( kind = 4 ) kk
integer ( kind = 4 ) lend(*)
integer ( kind = 4 ) list(*)
integer ( kind = 4 ) lnew
integer ( kind = 4 ) lp
integer ( kind = 4 ) lptr(*)
integer ( kind = 4 ) lsav
integer ( kind = 4 ) n0
integer ( kind = 4 ) next
integer ( kind = 4 ) nst
k = kk
nst = n0
!
! Traverse the boundary in clockwise order, inserting K as
! the first neighbor of each boundary node, and converting
! the boundary node to an interior node.
!
next = nst
do
lp = lend(next)
call insert ( k, lp, list, lptr, lnew )
next = -list(lp)
list(lp) = next
if ( next == nst ) then
exit
end if
end do
!
! Traverse the boundary again, adding each node to K's adjacency list.
!
lsav = lnew
do
lp = lend(next)
list(lnew) = next
lptr(lnew) = lnew + 1
lnew = lnew + 1
next = list(lp)
if ( next == nst ) then
exit
end if
end do
lptr(lnew-1) = lsav
lend(k) = lnew - 1
return
end
subroutine crlist ( n, ncol, x, y, z, list, lend, lptr, lnew, &
ltri, listc, nb, xc, yc, zc, rc, ier )
!*****************************************************************************80
!
!! CRLIST returns triangle circumcenters and other information.
!
! Discussion:
!
! Given a Delaunay triangulation of nodes on the surface
! of the unit sphere, this subroutine returns the set of
! triangle circumcenters corresponding to Voronoi vertices,
! along with the circumradii and a list of triangle indexes
! LISTC stored in one-to-one correspondence with LIST/LPTR
! entries.
!
! A triangle circumcenter is the point (unit vector) lying
! at the same angular distance from the three vertices and
! contained in the same hemisphere as the vertices. (Note
! that the negative of a circumcenter is also equidistant
! from the vertices.) If the triangulation covers the
! surface, the Voronoi vertices are the circumcenters of the
! triangles in the Delaunay triangulation. LPTR, LEND, and
! LNEW are not altered in this case.
!
! On the other hand, if the nodes are contained in a
! single hemisphere, the triangulation is implicitly extended
! to the entire surface by adding pseudo-arcs (of length
! greater than 180 degrees) between boundary nodes forming
! pseudo-triangles whose 'circumcenters' are included in the
! list. This extension to the triangulation actually
! consists of a triangulation of the set of boundary nodes in
! which the swap test is reversed (a non-empty circumcircle
! test). The negative circumcenters are stored as the
! pseudo-triangle 'circumcenters'. LISTC, LPTR, LEND, and
! LNEW contain a data structure corresponding to the
! extended triangulation (Voronoi diagram), but LIST is not
! altered in this case. Thus, if it is necessary to retain
! the original (unextended) triangulation data structure,
! copies of LPTR and LNEW must be saved before calling this
! routine.
!
! Modified:
!
! 16 June 2007
!
! Author:
!
! Robert Renka
!
! Reference:
!
! Robert Renka,
! Algorithm 772: STRIPACK,
! Delaunay Triangulation and Voronoi Diagram on the Surface of a Sphere,
! ACM Transactions on Mathematical Software,
! Volume 23, Number 3, September 1997, pages 416-434.
!
! Parameters:
!
! Input, integer ( kind = 4 ) N, the number of nodes in the triangulation.
! 3 <= N. Note that, if N = 3, there are only two Voronoi vertices
! separated by 180 degrees, and the Voronoi regions are not well defined.
!
! Input, integer ( kind = 4 ) NCOL, the number of columns reserved for LTRI.
! This must be at least NB-2, where NB is the number of boundary nodes.
!
! Input, real ( kind = 8 ) X(N), Y(N), Z(N), the coordinates of the nodes
! (unit vectors).
!
! Input, integer ( kind = 4 ) LIST(6*(N-2)), the set of adjacency lists.
! Refer to TRMESH.
!