-
Notifications
You must be signed in to change notification settings - Fork 35
/
surface_intersection.f90
1524 lines (1330 loc) · 60.7 KB
/
surface_intersection.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
! Licensed under the Apache License, Version 2.0 (the "License");
! you may not use this file except in compliance with the License.
! You may obtain a copy of the License at
!
! https://www.apache.org/licenses/LICENSE-2.0
!
! Unless required by applicable law or agreed to in writing, software
! distributed under the License is distributed on an "AS IS" BASIS,
! WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
! See the License for the specific language governing permissions and
! limitations under the License.
module surface_intersection
use, intrinsic :: iso_c_binding, only: c_double, c_int, c_bool
use status, only: &
Status_SUCCESS, Status_INSUFFICIENT_SPACE, Status_SAME_CURVATURE, &
Status_BAD_TANGENT, Status_EDGE_END, Status_UNKNOWN
use curve, only: CurveData, LOCATE_MISS, evaluate_hodograph, get_curvature
use curve_intersection, only: &
BoxIntersectionType_INTERSECTION, INTERSECTIONS_WORKSPACE, &
bbox_intersect, all_intersections
use helpers, only: cross_product, contains_nd, vector_close
use types, only: dp
use surface, only: &
evaluate_barycentric, jacobian_both, subdivide_nodes, compute_edge_nodes
implicit none
private &
LocateCandidate, MAX_LOCATE_SUBDIVISIONS, LOCATE_EPS, MAX_EDGES, &
newton_refine_solve, split_candidate, allocate_candidates, &
update_candidates, ignored_edge_corner, ignored_double_corner, &
ignored_corner, classify_tangent_intersection, no_intersections, &
remove_node, finalize_segment, check_contained
public &
Intersection, CurvedPolygonSegment, &
IntersectionClassification_FIRST, IntersectionClassification_SECOND, &
IntersectionClassification_OPPOSED, &
IntersectionClassification_TANGENT_FIRST, &
IntersectionClassification_TANGENT_SECOND, &
IntersectionClassification_IGNORED_CORNER, SurfaceContained_NEITHER, &
SurfaceContained_FIRST, SurfaceContained_SECOND, newton_refine, &
locate_point, classify_intersection, add_st_vals, &
surfaces_intersection_points, get_next, to_front, add_segment, &
interior_combine, surfaces_intersect, surfaces_intersect_abi, &
free_surface_intersections_workspace
! NOTE: This (for now) is not meant to be C-interoperable.
type :: Intersection
real(c_double) :: s = -1.0_dp
real(c_double) :: t = -1.0_dp
integer(c_int) :: index_first = -1
integer(c_int) :: index_second = -1
integer(c_int) :: interior_curve = -99 ! Hopefully an unused enum value
end type Intersection
! For ``locate_point``; this is not meant to be C-interoperable.
type :: LocateCandidate
real(c_double) :: centroid_x = 1.0_dp ! Actually triple.
real(c_double) :: centroid_y = 1.0_dp ! Actually triple.
real(c_double) :: width = 1.0_dp
real(c_double), allocatable :: nodes(:, :)
end type LocateCandidate
type, bind(c) :: CurvedPolygonSegment
real(c_double) :: start
real(c_double) :: end_
integer(c_int) :: edge_index
end type CurvedPolygonSegment
! NOTE: These values are also defined in equivalent Python source.
integer(c_int), parameter :: MAX_LOCATE_SUBDIVISIONS = 20
real(c_double), parameter :: LOCATE_EPS = 0.5_dp**47
integer(c_int), parameter :: MAX_EDGES = 10
! Values of IntersectionClassification enum:
integer(c_int), parameter :: IntersectionClassification_FIRST = 0
integer(c_int), parameter :: IntersectionClassification_SECOND = 1
integer(c_int), parameter :: IntersectionClassification_OPPOSED = 2
integer(c_int), parameter :: IntersectionClassification_TANGENT_FIRST = 3
integer(c_int), parameter :: IntersectionClassification_TANGENT_SECOND = 4
integer(c_int), parameter :: IntersectionClassification_IGNORED_CORNER = 5
! Values of SurfaceContained enum:
integer(c_int), parameter :: SurfaceContained_NEITHER = 0
integer(c_int), parameter :: SurfaceContained_FIRST = 1
integer(c_int), parameter :: SurfaceContained_SECOND = 2
! Long-lived workspaces for ``surfaces_intersect_abi()``. If multiple
! threads are used, each of these **should** be thread-local.
integer(c_int), allocatable :: SEGMENT_ENDS_WORKSPACE(:)
type(CurvedPolygonSegment), allocatable :: SEGMENTS_WORKSPACE(:)
contains
subroutine newton_refine_solve( &
jac_both, x_val, surf_x, y_val, surf_y, delta_s, delta_t)
! NOTE: This is a helper for ``newton_refine``.
real(c_double), intent(in) :: jac_both(1, 4)
real(c_double), intent(in) :: x_val, surf_x
real(c_double), intent(in) :: y_val, surf_y
real(c_double), intent(out) :: delta_s, delta_t
! Variables outside of signature.
real(c_double) :: e_val, f_val, denominator
e_val = x_val - surf_x
f_val = y_val - surf_y
denominator = ( &
jac_both(1, 1) * jac_both(1, 4) - jac_both(1, 2) * jac_both(1, 3))
delta_s = (jac_both(1, 4) * e_val - jac_both(1, 3) * f_val) / denominator
delta_t = (jac_both(1, 1) * f_val - jac_both(1, 2) * e_val) / denominator
end subroutine newton_refine_solve
subroutine newton_refine( &
num_nodes, nodes, degree, x_val, y_val, &
s, t, updated_s, updated_t) &
bind(c, name='newton_refine_surface')
integer(c_int), intent(in) :: num_nodes
real(c_double), intent(in) :: nodes(num_nodes, 2)
integer(c_int), intent(in) :: degree
real(c_double), intent(in) :: x_val, y_val
real(c_double), intent(in) :: s, t
real(c_double), intent(out) :: updated_s, updated_t
! Variables outside of signature.
real(c_double) :: point(1, 2), jac_both(1, 4)
real(c_double) :: jac_nodes(num_nodes - degree - 1, 4)
real(c_double) :: delta_s, delta_t
call evaluate_barycentric( &
num_nodes, 2, nodes, degree, &
1.0_dp - s - t, s, t, point)
if (point(1, 1) == x_val .AND. point(1, 2) == y_val) then
! No refinement is needed.
updated_s = s
updated_t = t
return
end if
call jacobian_both( &
num_nodes, 2, nodes, degree, jac_nodes)
call evaluate_barycentric( &
num_nodes - degree - 1, 4, jac_nodes, degree - 1, &
1.0_dp - s - t, s, t, jac_both)
call newton_refine_solve( &
jac_both, x_val, point(1, 1), &
y_val, point(1, 2), delta_s, delta_t)
updated_s = s + delta_s
updated_t = t + delta_t
end subroutine newton_refine
subroutine split_candidate( &
num_nodes, degree, candidate, num_next_candidates, next_candidates)
! NOTE: This assumes nodes are 2-dimensional.
! NOTE: This assumes that the nodes in each sub-candidate are
! not yet allocated.
integer(c_int), intent(in) :: num_nodes
integer(c_int), intent(in) :: degree
type(LocateCandidate), intent(in) :: candidate
integer(c_int), intent(in) :: num_next_candidates
type(LocateCandidate), intent(inout) :: next_candidates(:)
! Variables outside of signature.
real(c_double) :: half_width
! Allocate the new nodes and call sub-divide.
! NOTE: This **assumes** but does not check that if the nodes are
! allocated, they are also the correct shape.
if (.NOT. allocated(next_candidates(num_next_candidates - 3)%nodes)) then
allocate(next_candidates(num_next_candidates - 3)%nodes(num_nodes, 2))
end if
if (.NOT. allocated(next_candidates(num_next_candidates - 2)%nodes)) then
allocate(next_candidates(num_next_candidates - 2)%nodes(num_nodes, 2))
end if
if (.NOT. allocated(next_candidates(num_next_candidates - 1)%nodes)) then
allocate(next_candidates(num_next_candidates - 1)%nodes(num_nodes, 2))
end if
if (.NOT. allocated(next_candidates(num_next_candidates)%nodes)) then
allocate(next_candidates(num_next_candidates)%nodes(num_nodes, 2))
end if
call subdivide_nodes( &
num_nodes, 2, candidate%nodes, degree, &
next_candidates(num_next_candidates - 3)%nodes, &
next_candidates(num_next_candidates - 2)%nodes, &
next_candidates(num_next_candidates - 1)%nodes, &
next_candidates(num_next_candidates)%nodes)
half_width = 0.5_dp * candidate%width
! Subdivision A.
next_candidates(num_next_candidates - 3)%centroid_x = ( &
candidate%centroid_x - half_width)
next_candidates(num_next_candidates - 3)%centroid_y = ( &
candidate%centroid_y - half_width)
next_candidates(num_next_candidates - 3)%width = half_width
! Subdivision B.
next_candidates(num_next_candidates - 2)%centroid_x = candidate%centroid_x
next_candidates(num_next_candidates - 2)%centroid_y = candidate%centroid_y
next_candidates(num_next_candidates - 2)%width = -half_width
! Subdivision C.
next_candidates(num_next_candidates - 1)%centroid_x = ( &
candidate%centroid_x + candidate%width)
next_candidates(num_next_candidates - 1)%centroid_y = ( &
next_candidates(num_next_candidates - 3)%centroid_y)
next_candidates(num_next_candidates - 1)%width = half_width
! Subdivision D.
next_candidates(num_next_candidates)%centroid_x = ( &
next_candidates(num_next_candidates - 3)%centroid_x)
next_candidates(num_next_candidates)%centroid_y = ( &
candidate%centroid_y + candidate%width)
next_candidates(num_next_candidates)%width = half_width
end subroutine split_candidate
subroutine allocate_candidates(num_candidates, next_candidates)
integer(c_int), intent(in) :: num_candidates
type(LocateCandidate), allocatable, intent(inout) :: next_candidates(:)
if (allocated(next_candidates)) then
if (size(next_candidates) < 4 * num_candidates) then
deallocate(next_candidates)
allocate(next_candidates(4 * num_candidates))
end if
else
allocate(next_candidates(4 * num_candidates))
end if
end subroutine allocate_candidates
subroutine update_candidates( &
num_nodes, degree, point, num_candidates, candidates, &
num_next_candidates, next_candidates)
integer(c_int), intent(in) :: num_nodes, degree
real(c_double), intent(in) :: point(1, 2)
integer(c_int), intent(in) :: num_candidates
type(LocateCandidate), intent(in) :: candidates(:)
integer(c_int), intent(inout) :: num_next_candidates
type(LocateCandidate), allocatable, intent(inout) :: next_candidates(:)
! Variables outside of signature.
integer(c_int) :: cand_index
logical(c_bool) :: predicate
call allocate_candidates(num_candidates, next_candidates)
do cand_index = 1, num_candidates
call contains_nd( &
num_nodes, 2, candidates(cand_index)%nodes, &
point(1, :), predicate)
if (predicate) then
num_next_candidates = num_next_candidates + 4
call split_candidate( &
num_nodes, degree, candidates(cand_index), &
num_next_candidates, next_candidates)
end if
end do
end subroutine update_candidates
subroutine locate_point( &
num_nodes, nodes, degree, x_val, y_val, s_val, t_val) &
bind(c, name='locate_point_surface')
! NOTE: This solves the inverse problem B(s, t) = (x, y) (if it can be
! solved). Does so by subdividing the surface until the sub-surfaces
! are sufficiently small, then using Newton's method to narrow
! in on the pre-image of the point.
! NOTE: This returns ``-1`` (``LOCATE_MISS``) for ``s_val`` as a signal
! for "point is not on the surface".
! NOTE: This assumes, but does not check, that the surface is "valid",
! i.e. all pre-image are unique.
integer(c_int), intent(in) :: num_nodes
real(c_double), intent(in) :: nodes(num_nodes, 2)
integer(c_int), intent(in) :: degree
real(c_double), intent(in) :: x_val, y_val
real(c_double), intent(out) :: s_val, t_val
! Variables outside of signature.
real(c_double) :: point(1, 2)
integer(c_int) :: sub_index
integer(c_int) :: num_candidates, num_next_candidates
type(LocateCandidate), allocatable :: candidates_odd(:), candidates_even(:)
real(c_double) :: s_approx, t_approx
real(c_double) :: actual(1, 2)
logical(c_bool) :: is_even
point(1, :) = [x_val, y_val]
! Start out with the full curve.
allocate(candidates_odd(1)) ! First iteration is odd.
candidates_odd(1) = LocateCandidate(1.0_dp, 1.0_dp, 1.0_dp, nodes)
num_candidates = 1
s_val = LOCATE_MISS
is_even = .TRUE. ! At zero.
do sub_index = 1, MAX_LOCATE_SUBDIVISIONS + 1
is_even = .NOT. is_even ! Switch parity.
num_next_candidates = 0
if (is_even) then
! Read from even, write to odd.
call update_candidates( &
num_nodes, degree, point, num_candidates, candidates_even, &
num_next_candidates, candidates_odd)
else
! Read from odd, write to even.
call update_candidates( &
num_nodes, degree, point, num_candidates, candidates_odd, &
num_next_candidates, candidates_even)
end if
! If there are no more candidates, we are done.
if (num_next_candidates == 0) then
return
end if
num_candidates = num_next_candidates
end do
! Compute the s- and t-parameters as the mean of the centroid positions.
! We divide by `3n` rather than `n` since we have tracked thrice the
! centroid rather than the centroid itself (to avoid round-off until
! right now).
if (is_even) then
! NOTE: We "exclude" this block from ``lcov`` because it can never
! be invoked while ``MAX_LOCATE_SUBDIVISIONS`` is even.
! LCOV_EXCL_START
s_approx = ( &
sum(candidates_odd(:num_candidates)%centroid_x) / &
(3.0_dp * num_candidates))
t_approx = ( &
sum(candidates_odd(:num_candidates)%centroid_y) / &
(3.0_dp * num_candidates))
! LCOV_EXCL_STOP
else
s_approx = ( &
sum(candidates_even(:num_candidates)%centroid_x) / &
(3.0_dp * num_candidates))
t_approx = ( &
sum(candidates_even(:num_candidates)%centroid_y) / &
(3.0_dp * num_candidates))
end if
call newton_refine( &
num_nodes, nodes, degree, x_val, y_val, &
s_approx, t_approx, s_val, t_val)
! Check if the solution is "close enough" ...
call evaluate_barycentric( &
num_nodes, 2, nodes, degree, &
1.0_dp - s_val - t_val, s_val, t_val, actual)
if (vector_close(2, point, actual, LOCATE_EPS)) then
return
end if
! ... and if **not** close enough, do one more round of
! Newton's method if not.
s_approx = s_val
t_approx = t_val
call newton_refine( &
num_nodes, nodes, degree, x_val, y_val, &
s_approx, t_approx, s_val, t_val)
end subroutine locate_point
logical(c_bool) function ignored_edge_corner( &
edge_tangent, corner_tangent, num_nodes, &
previous_edge_nodes) result(predicate)
real(c_double), intent(in) :: edge_tangent(1, 2)
real(c_double), intent(in) :: corner_tangent(1, 2)
integer(c_int), intent(in) :: num_nodes
real(c_double), intent(in) :: previous_edge_nodes(num_nodes, 2)
! Variables outside of signature.
real(c_double) :: cross_prod
real(c_double) :: alt_corner_tangent(1, 2)
call cross_product( &
edge_tangent, corner_tangent, cross_prod)
! A negative cross product indicates that ``edge_tangent`` is
! "inside" / "to the left" of ``corner_tangent`` (due to right-hand rule).
if (cross_prod > 0.0_dp) then
predicate = .FALSE.
return
end if
! Do the same for the **other** tangent at the corner.
call evaluate_hodograph( &
1.0_dp, num_nodes, 2, &
previous_edge_nodes, alt_corner_tangent)
! Change the direction of the "in" tangent so that it points "out".
alt_corner_tangent = -alt_corner_tangent
call cross_product( &
edge_tangent, alt_corner_tangent, cross_prod)
predicate = (cross_prod <= 0.0_dp)
end function ignored_edge_corner
logical(c_bool) function ignored_double_corner( &
edges_first, edges_second, intersection_, &
tangent_s, tangent_t) result(predicate)
! NOTE: This assumes that ``intersection_%index_(first|second)``
! are in [1, 2, 3] (marked as the edges of a surface).
type(CurveData), intent(in) :: edges_first(3), edges_second(3)
type(Intersection), intent(in) :: intersection_
real(c_double), intent(in) :: tangent_s(1, 2), tangent_t(1, 2)
! Variables outside of signature.
integer(c_int) :: index, num_nodes
real(c_double) :: alt_tangent_s(1, 2), alt_tangent_t(1, 2)
real(c_double) :: cross_prod1, cross_prod2, cross_prod3
! Compute the other edge for the ``s`` surface.
index = 1 + modulo(intersection_%index_first - 2, 3)
num_nodes = size(edges_first(index)%nodes, 1)
call evaluate_hodograph( &
1.0_dp, num_nodes, 2, &
edges_first(index)%nodes, alt_tangent_s)
! First check if ``tangent_t`` is interior to the ``s`` surface.
call cross_product( &
tangent_s, tangent_t, cross_prod1)
! A positive cross product indicates that ``tangent_t`` is
! interior to ``tangent_s``. Similar for ``alt_tangent_s``.
! If ``tangent_t`` is interior to both, then the surfaces
! do more than just "kiss" at the corner, so the corner should
! not be ignored.
if (cross_prod1 >= 0.0_dp) then
! Only compute ``cross_prod2`` if we need to.
call cross_product( &
alt_tangent_s, tangent_t, cross_prod2)
if (cross_prod2 >= 0.0_dp) then
predicate = .FALSE.
return
end if
end if
! If ``tangent_t`` is not interior, we check the other ``t``
! edge that ends at the corner.
index = 1 + modulo(intersection_%index_second - 2, 3)
num_nodes = size(edges_second(index)%nodes, 1)
call evaluate_hodograph( &
1.0_dp, num_nodes, 2, &
edges_second(index)%nodes, alt_tangent_t)
! Change the direction of the "in" tangent so that it points "out".
alt_tangent_t = -alt_tangent_t
call cross_product( &
tangent_s, alt_tangent_t, cross_prod2)
if (cross_prod2 >= 0.0_dp) then
! Only compute ``cross_prod3`` if we need to.
call cross_product( &
alt_tangent_s, alt_tangent_t, cross_prod3)
if (cross_prod3 >= 0.0_dp) then
predicate = .FALSE.
return
end if
end if
! If neither of ``tangent_t`` or ``alt_tangent_t`` are interior
! to the ``s`` surface, one of two things is true. Either
! the two surfaces have no interior intersection (1) or the
! ``s`` surface is bounded by both edges of the ``t`` surface
! at the corner intersection (2). To detect (2), we only need
! check if ``tangent_s`` is interior to both ``tangent_t``
! and ``alt_tangent_t``. ``cross_prod1`` contains
! (tangent_s) x (tangent_t), so it's negative will tell if
! ``tangent_s`` is interior. Similarly, ``cross_prod3``
! contains (tangent_s) x (alt_tangent_t), but we also reversed
! the sign on ``alt_tangent_t`` so switching the sign back
! and reversing the arguments in the cross product cancel out.
predicate = (cross_prod1 > 0.0_dp .OR. cross_prod2 < 0.0_dp)
end function ignored_double_corner
logical(c_bool) function ignored_corner( &
edges_first, edges_second, intersection_, &
tangent_s, tangent_t) result(predicate)
! NOTE: This assumes that ``intersection_%index_(first|second)``
! are in [1, 2, 3] (marked as the edges of a surface).
type(CurveData), intent(in) :: edges_first(3), edges_second(3)
type(Intersection), intent(in) :: intersection_
real(c_double), intent(in) :: tangent_s(1, 2), tangent_t(1, 2)
! Variables outside of signature.
integer(c_int) :: index, num_nodes
if (intersection_%s == 0.0_dp) then
if (intersection_%t == 0.0_dp) then
! Double corner.
predicate = ignored_double_corner( &
edges_first, edges_second, intersection_, &
tangent_s, tangent_t)
else
! s-only corner.
index = 1 + modulo(intersection_%index_first - 2, 3)
! NOTE: This is a "trick" which requires `modulo` (not `mod`)
! since it will return values in {0, 1, 2}. The goal is
! to send [1, 2, 3] --> [3, 1, 2] (i.e. move indices to
! the left and wrap around).
num_nodes = size(edges_first(index)%nodes, 1)
predicate = ignored_edge_corner( &
tangent_t, tangent_s, num_nodes, edges_first(index)%nodes)
end if
else if (intersection_%t == 0.0_dp) then
! t-only corner.
index = 1 + modulo(intersection_%index_second - 2, 3)
! NOTE: This is a "trick" which requires `modulo` (not `mod`)
! since it will return values in {0, 1, 2}. The goal is
! to send [1, 2, 3] --> [3, 1, 2] (i.e. move indices to
! the left and wrap around).
num_nodes = size(edges_second(index)%nodes, 1)
predicate = ignored_edge_corner( &
tangent_s, tangent_t, num_nodes, edges_second(index)%nodes)
else
predicate = .FALSE.
end if
end function ignored_corner
subroutine classify_tangent_intersection( &
edges_first, edges_second, intersection_, &
tangent_s, tangent_t, enum_, status)
! NOTE: This **assumes**, but does not check that
! ``intersection_%index_(first|second)`` are valid indices within
! ``edges_(first|second)`` and that each of those ``CurveData``
! instances have already allocated ``%nodes``.
! Possible error states:
! * Status_SUCCESS : On success.
! * Status_BAD_TANGENT : If the curves move in an opposite direction
! while defining overlapping arcs.
! * Status_SAME_CURVATURE: If the curves have identical curvature at the
! point of tangency.
type(CurveData), intent(in) :: edges_first(3), edges_second(3)
type(Intersection), intent(in) :: intersection_
real(c_double), intent(in) :: tangent_s(1, 2), tangent_t(1, 2)
integer(c_int), intent(out) :: enum_
integer(c_int), intent(out) :: status
! Variables outside of signature.
real(c_double) :: dot_prod
integer(c_int) :: num_nodes
real(c_double) :: curvature1, curvature2
real(c_double) :: sign1, sign2
real(c_double) :: delta_c
status = Status_SUCCESS
dot_prod = dot_product(tangent_s(1, :), tangent_t(1, :))
! NOTE: When computing curvatures we assume that we don't have lines
! here, because lines that are tangent at an intersection are
! parallel and we don't handle that case.
num_nodes = size(edges_first(intersection_%index_first)%nodes, 1)
call get_curvature( &
num_nodes, 2, edges_first(intersection_%index_first)%nodes, &
tangent_s, intersection_%s, curvature1)
num_nodes = size(edges_second(intersection_%index_second)%nodes, 1)
call get_curvature( &
num_nodes, 2, edges_second(intersection_%index_second)%nodes, &
tangent_t, intersection_%t, curvature2)
if (dot_prod < 0.0_dp) then
! If the tangent vectors are pointing in the opposite direction,
! then the curves are facing opposite directions.
sign1 = sign(1.0_dp, curvature1)
sign2 = sign(1.0_dp, curvature2)
if (sign1 == sign2) then
! If both curvatures are positive, since the curves are
! moving in opposite directions, the tangency isn't part of
! the surface intersection.
if (sign1 == 1.0_dp) then
enum_ = IntersectionClassification_OPPOSED
else
status = Status_BAD_TANGENT
end if
else
delta_c = abs(curvature1) - abs(curvature2)
if (delta_c == 0.0_dp) then
status = Status_SAME_CURVATURE
else
sign2 = sign(1.0_dp, delta_c)
if (sign1 == sign2) then
enum_ = IntersectionClassification_OPPOSED
else
status = Status_BAD_TANGENT
end if
end if
end if
else
if (curvature1 > curvature2) then
enum_ = IntersectionClassification_TANGENT_FIRST
else if (curvature1 < curvature2) then
enum_ = IntersectionClassification_TANGENT_SECOND
else
status = Status_SAME_CURVATURE
end if
end if
end subroutine classify_tangent_intersection
subroutine classify_intersection( &
edges_first, edges_second, intersection_, enum_, status)
! NOTE: This is **explicitly** not intended for C inter-op.
! NOTE: This subroutine is not part of the C ABI for this module,
! but it is (for now) public, so that it can be tested.
! NOTE: This **assumes**, but does not check that
! ``intersection_%index_(first|second)`` are valid indices within
! ``edges_(first|second)`` and that each of those ``CurveData``
! instances have already allocated ``%nodes``.
! Possible error states:
! * Status_SUCCESS : On success.
! * Status_EDGE_END : If either the s- or t-parameter of the
! intersection is equal to ``1.0``, i.e. if the
! intersection is at the end of an edge. Since
! another edge of the surface will also contain
! that point (at it's beginning), that edge
! should be used.
! * Status_BAD_TANGENT : Via ``classify_tangent_intersection()``.
! * Status_SAME_CURVATURE: Via ``classify_tangent_intersection()``.
type(CurveData), intent(in) :: edges_first(3), edges_second(3)
type(Intersection), intent(in) :: intersection_
integer(c_int), intent(out) :: enum_
integer(c_int), intent(out) :: status
! Variables outside of signature.
integer(c_int) :: num_nodes
real(c_double) :: tangent_s(1, 2), tangent_t(1, 2)
real(c_double) :: cross_prod
status = Status_SUCCESS
if (intersection_%s == 1.0_dp .OR. intersection_%t == 1.0_dp) then
status = Status_EDGE_END
return
end if
! NOTE: We assume, but don't check that ``%nodes`` is allocated
! and that it has 2 columns.
num_nodes = size(edges_first(intersection_%index_first)%nodes, 1)
call evaluate_hodograph( &
intersection_%s, num_nodes, 2, &
edges_first(intersection_%index_first)%nodes, tangent_s)
! NOTE: We assume, but don't check that ``%nodes`` is allocated
! and that it has 2 columns.
num_nodes = size(edges_second(intersection_%index_second)%nodes, 1)
call evaluate_hodograph( &
intersection_%t, num_nodes, 2, &
edges_second(intersection_%index_second)%nodes, tangent_t)
if (ignored_corner( &
edges_first, edges_second, intersection_, &
tangent_s, tangent_t)) then
enum_ = IntersectionClassification_IGNORED_CORNER
return
end if
! Take the cross product of tangent vectors to determine which one
! is more "inside" / "to the left".
call cross_product( &
tangent_s, tangent_t, cross_prod)
if (cross_prod < 0.0_dp) then
enum_ = IntersectionClassification_FIRST
else if (cross_prod > 0.0_dp) then
enum_ = IntersectionClassification_SECOND
else
call classify_tangent_intersection( &
edges_first, edges_second, intersection_, &
tangent_s, tangent_t, enum_, status)
end if
end subroutine classify_intersection
subroutine add_st_vals( &
edges_first, edges_second, num_st_vals, st_vals, &
index_first, index_second, intersections, num_intersections, &
all_types, status)
! NOTE: This is **explicitly** not intended for C inter-op.
! NOTE: This subroutine is not part of the C ABI for this module,
! but it is (for now) public, so that it can be tested.
! NOTE: This assumes but does not check that ``num_st_vals > 0``.
! NOTE: This assumes but does not check that each of ``index_first``
! and ``index_second`` is in {1, 2, 3}.
! Possible error states:
! * Status_SUCCESS : On success.
! * Status_EDGE_END : Via ``classify_intersection()``.
! * Status_BAD_TANGENT : Via ``classify_intersection()``.
! * Status_SAME_CURVATURE: Via ``classify_intersection()``.
type(CurveData), intent(in) :: edges_first(3), edges_second(3)
integer(c_int), intent(in) :: num_st_vals
real(c_double), intent(in) :: st_vals(2, num_st_vals)
integer(c_int), intent(in) :: index_first, index_second
type(Intersection), allocatable, intent(inout) :: intersections(:)
integer(c_int), intent(inout) :: num_intersections
integer(c_int), intent(inout) :: all_types
integer(c_int), intent(out) :: status
! Variables outside of signature.
integer(c_int) :: curr_size, i, intersection_index
type(Intersection), allocatable :: intersections_swap(:)
integer(c_int) :: enum_
status = Status_SUCCESS
intersection_index = num_intersections
! NOTE: Intersections at the end of an edge will be skipped, so this
! may over-estimate the number of intersections.
num_intersections = num_intersections + num_st_vals
if (allocated(intersections)) then
curr_size = size(intersections)
if (curr_size < num_intersections) then
allocate(intersections_swap(num_intersections))
intersections_swap(:curr_size) = intersections(:curr_size)
call move_alloc(intersections_swap, intersections)
end if
else
allocate(intersections(num_intersections))
end if
do i = 1, num_st_vals
if (st_vals(1, i) == 1.0_dp .OR. st_vals(2, i) == 1.0_dp) then
! If the intersection is at the end of an edge, ignore it.
! This intersection **could** be saved to verify that it matches an
! equivalent intersection from the beginning of the previous edge.
cycle
end if
intersection_index = intersection_index + 1
intersections(intersection_index)%s = st_vals(1, i)
intersections(intersection_index)%t = st_vals(2, i)
intersections(intersection_index)%index_first = index_first
intersections(intersection_index)%index_second = index_second
call classify_intersection( &
edges_first, edges_second, &
intersections(intersection_index), enum_, status)
if (status /= Status_SUCCESS) then
return
end if
! NOTE: This assumes, but does not check, that ``enum_`` is in
! {0, 1, 2, 3, 4, 5} (should be an enum value from
! ``IntersectionClassification``) which limits the value of
! ``all_types`` to [0, 63] (inclusive).
all_types = ior(all_types, 2**enum_)
if ( &
enum_ == IntersectionClassification_FIRST .OR. &
enum_ == IntersectionClassification_SECOND) then
! "Keep" the intersection if it is ``FIRST`` or ``SECOND``.
intersections(intersection_index)%interior_curve = enum_
else
! "Discard" the intersection (rather, allow it to be over-written).
intersection_index = intersection_index - 1
end if
end do
! Actually set ``num_intersections`` based on the number of **accepted**
! ``s-t`` pairs.
num_intersections = intersection_index
end subroutine add_st_vals
subroutine surfaces_intersection_points( &
num_nodes1, nodes1, degree1, &
num_nodes2, nodes2, degree2, &
intersections, num_intersections, all_types, status)
! NOTE: This is **explicitly** not intended for C inter-op.
! NOTE: This (and ``add_st_vals``) ignores duplicate nodes (caused when
! an intersection happens at a corner of the surface, which will
! be at the end of one edge and the start of another). If desired,
! a "verify" option could be added (as is done in Python) to make
! sure that any "duplicate" intersections (i.e. ones that occur at
! the end of an edge) do match with the corresponding intersection
! at the beginning of an edge (in the case of a double corner, this
! means **three** duplicate intersections).
! NOTE: The returned ``all_types`` uses the first 6 bits to identify
! which classified states are among the intersections. (The
! possible classifications must be in {0, 1, 2, 3, 4, 5}).
! Possible error states:
! * Status_SUCCESS : On success.
! * Status_PARALLEL : Via ``curve_intersection.all_intersections()``.
! * Status_WIGGLE_FAIL : Via ``curve_intersection.all_intersections()``.
! * Status_NO_CONVERGE : Via ``curve_intersection.all_intersections()``.
! * (N >= MAX_CANDIDATES): Via ``curve_intersection.all_intersections()``.
! * Status_EDGE_END : Via ``add_st_vals()``.
! * Status_BAD_TANGENT : Via ``add_st_vals()``.
! * Status_SAME_CURVATURE: Via ``add_st_vals()``.
integer(c_int), intent(in) :: num_nodes1
real(c_double), intent(in) :: nodes1(num_nodes1, 2)
integer(c_int), intent(in) :: degree1
integer(c_int), intent(in) :: num_nodes2
real(c_double), intent(in) :: nodes2(num_nodes2, 2)
integer(c_int), intent(in) :: degree2
type(Intersection), allocatable, intent(inout) :: intersections(:)
integer(c_int), intent(out) :: num_intersections
integer(c_int), intent(out) :: all_types
integer(c_int), intent(out) :: status
! Variables outside of signature.
type(CurveData) :: edges_first(3), edges_second(3)
integer(c_int) :: index1, index2
integer(c_int) :: num_st_vals
! Compute the edge nodes for the first surface.
allocate(edges_first(1)%nodes(degree1 + 1, 2))
allocate(edges_first(2)%nodes(degree1 + 1, 2))
allocate(edges_first(3)%nodes(degree1 + 1, 2))
call compute_edge_nodes( &
num_nodes1, 2, nodes1, degree1, &
edges_first(1)%nodes, edges_first(2)%nodes, edges_first(3)%nodes)
! Compute the edge nodes for the second surface.
allocate(edges_second(1)%nodes(degree2 + 1, 2))
allocate(edges_second(2)%nodes(degree2 + 1, 2))
allocate(edges_second(3)%nodes(degree2 + 1, 2))
call compute_edge_nodes( &
num_nodes2, 2, nodes2, degree2, &
edges_second(1)%nodes, edges_second(2)%nodes, edges_second(3)%nodes)
num_intersections = 0
all_types = 0
do index1 = 1, 3
do index2 = 1, 3
call all_intersections( &
degree1 + 1, edges_first(index1)%nodes, &
degree2 + 1, edges_second(index2)%nodes, &
INTERSECTIONS_WORKSPACE, num_st_vals, status)
if (status == Status_SUCCESS) then
if (num_st_vals > 0) then
call add_st_vals( &
edges_first, edges_second, &
num_st_vals, INTERSECTIONS_WORKSPACE(:, :num_st_vals), &
index1, index2, intersections, &
num_intersections, all_types, status)
if (status /= Status_SUCCESS) then
return
end if
end if
else
return
end if
end do
end do
end subroutine surfaces_intersection_points
subroutine no_intersections( &
num_nodes1, nodes1, degree1, &
num_nodes2, nodes2, degree2, contained)
! NOTE: This is a helper for ``surfaces_intersect()``.
integer(c_int), intent(in) :: num_nodes1
real(c_double), intent(in) :: nodes1(num_nodes1, 2)
integer(c_int), intent(in) :: degree1
integer(c_int), intent(in) :: num_nodes2
real(c_double), intent(in) :: nodes2(num_nodes2, 2)
integer(c_int), intent(in) :: degree2
integer(c_int), intent(out) :: contained
! Variables outside of signature.
real(c_double) :: s_val, t_val
! If the first corner of ``surface1`` is contained in ``surface2``,
! then the whole surface must be since there are no intersections.
call locate_point( &
num_nodes2, nodes2, degree2, &
nodes1(1, 1), nodes1(1, 2), s_val, t_val)
if (s_val /= LOCATE_MISS) then
contained = SurfaceContained_FIRST
return
end if
! If the first corner of ``surface2`` is contained in ``surface1``,
! then the whole surface must be since there are no intersections.
call locate_point( &
num_nodes1, nodes1, degree1, &
nodes2(1, 1), nodes2(1, 2), s_val, t_val)
if (s_val /= LOCATE_MISS) then
contained = SurfaceContained_SECOND
return
end if
! If neither corner is contained, then there is no intersection.
contained = SurfaceContained_NEITHER
end subroutine no_intersections
subroutine remove_node(node, values, remaining)
! Removes a value from a list, if it is contained there.
! For example, if ``node == 4``, ``values == [1, 2, 4, 7, ...]`` and
! ``remainining == 4``, then after this subroutine will update to
! ``values == [1, 2, 7, ...]`` and ``remainining == 3``.
! NOTE: This assumes the values in ``values`` are in ascending order
! and unique.
! NOTE: This assumes that ``remaining`` does not exceed ``size(values)``
integer(c_int), intent(in) :: node
integer(c_int), intent(inout) :: values(:)
integer(c_int), intent(inout) :: remaining
! Variables outside of signature.
integer(c_int) :: i
! NOTE: Since we assume the values are in ascending order and unique, a
! binary search could speed this up. In practice, ``values``
! will likely be too small for this to matter.
do i = 1, remaining
if (values(i) == node) then
! Remove ``values(i)`` by shifting down all remaining values.
values(i:remaining - 1) = values(i + 1:remaining)
remaining = remaining - 1
return
end if
end do
end subroutine remove_node
subroutine get_next( &
num_intersections, intersections, unused, remaining, &
start, curr_node, next_node, at_start)
! NOTE: This subroutine is not meant to be part of the interface for this
! module, but it is (for now) public, so that it can be tested.
integer(c_int), intent(in) :: num_intersections
type(Intersection), intent(in) :: intersections(num_intersections)
! ``unused`` contains the indices of intersections that have not yet been
! used as a node.
integer(c_int), intent(inout) :: unused(:)
integer(c_int), intent(inout) :: remaining
integer(c_int), intent(in) :: start
type(Intersection), intent(in) :: curr_node
type(Intersection), intent(out) :: next_node
logical(c_bool), intent(out) :: at_start
! Variables outside of signature.
real(c_double) :: edge_param
integer(c_int) :: i, intersection_index
intersection_index = -1
at_start = .FALSE.
if (curr_node%interior_curve == IntersectionClassification_FIRST) then
do i = 1, num_intersections
if ( &
intersections(i)%index_first == curr_node%index_first .AND. &
intersections(i)%s > curr_node%s) then
if (intersection_index == -1) then
intersection_index = i
edge_param = intersections(i)%s
else
if (intersections(i)%s < edge_param) then
intersection_index = i
edge_param = intersections(i)%s
end if
end if
end if
end do
! If there is no other intersection on the edge, just return
! the segment end.
if (intersection_index == -1) then
! NOTE: We **explicitly** do not set ``t`` or ``index_second``.
next_node%index_first = curr_node%index_first
next_node%s = 1.0_dp
next_node%interior_curve = IntersectionClassification_FIRST
else
next_node = intersections(intersection_index)
at_start = (intersection_index == start)
! Remove the index from the set of ``unused`` intersections, if
! it is contained there.
call remove_node(intersection_index, unused, remaining)
end if
else
! NOTE: This assumes but does not check that
! ``curr_node%interior_curve`` is ``SECOND``.
do i = 1, num_intersections
if ( &
intersections(i)%index_second == curr_node%index_second .AND. &
intersections(i)%t > curr_node%t) then
if (intersection_index == -1) then
intersection_index = i
edge_param = intersections(i)%t
else
if (intersections(i)%t < edge_param) then
intersection_index = i
edge_param = intersections(i)%t
end if