forked from LadaF/Fortran---CGAL-polyhedra
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathcgal_polyhedra.f90
177 lines (142 loc) · 4.94 KB
/
cgal_polyhedra.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
!(C) Vladimir Fuka 2014
!GNU GPL license v3
module CGAL_Polyhedra
use iso_c_binding
implicit none
private
public cgal_polyhedron_read, &
cgal_polyhedron_closest, &
cgal_polyhedron_intersects_ray, &
cgal_polyhedron_bbox, &
cgal_polyhedron_finalize
type, bind(C) :: d3
real(c_double) :: x, y, z
end type
interface
subroutine polyhedron_from_file(ptree, fname, verbose, ierr) bind(C,name="polyhedron_from_file")
import
type(c_ptr),intent(out) :: ptree
character(kind=c_char),intent(in) :: fname(*)
integer(c_int),value :: verbose !avoid bool in C++, not equal to c_bool
integer(c_int),intent(out) :: ierr
end subroutine
subroutine polyhedron_closest(ptree, query, near) bind(C,name="polyhedron_closest")
import
type(c_ptr),value :: ptree
type(d3),intent(in) :: query
type(d3),intent(out) :: near
end subroutine
function polyhedron_intersects_ray(ptree, origin, vec) result(res) bind(C,name="polyhedron_intersects_ray")
import
logical(c_bool) :: res
type(c_ptr),value :: ptree
type(d3),intent(in) :: origin, vec
end function
subroutine polyhedron_bbox(ptree, min, max) bind(C,name="polyhedron_bbox")
import
type(c_ptr),value :: ptree
type(d3),intent(out) :: min, max
end subroutine
subroutine polyhedron_finalize(ptree) bind(C,name="polyhedron_finalize")
import
type(c_ptr),intent(inout) :: ptree
end subroutine
end interface
interface cgal_polyhedron_closest
module procedure cgal_polyhedron_closest_s
module procedure cgal_polyhedron_closest_d
end interface
interface cgal_polyhedron_intersects_ray
module procedure cgal_polyhedron_intersects_ray_s
module procedure cgal_polyhedron_intersects_ray_d
end interface
interface cgal_polyhedron_bbox
module procedure cgal_polyhedron_bbox_s
module procedure cgal_polyhedron_bbox_d
end interface
contains
subroutine cgal_polyhedron_read(ptree, fname)
type(c_ptr),intent(out) :: ptree
character(*),intent(in) :: fname
integer(c_int) :: ierr
integer :: imaster
call polyhedron_from_file(ptree, fname//c_null_char, 1, ierr)
if (ierr==1) then
write (*,*) "Error reading file "//fname//", it appears to be empty."
stop
else if (ierr==2) then
write (*,*) "Error reading file "//fname//"."
stop
end if
end subroutine
subroutine cgal_polyhedron_closest_s(ptree, xq,yq,zq, xn,yn,zn)
type(c_ptr),intent(in) :: ptree
real(c_float),intent(in) :: xq,yq,zq
real(c_float),intent(out) :: xn,yn,zn
type(d3) :: query, near
query = d3(real(xq,c_double),real(yq,c_double),real(zq,c_double))
call polyhedron_closest(ptree, query, near)
xn = real(near%x,c_float)
yn = real(near%y,c_float)
zn = real(near%z,c_float)
end subroutine
subroutine cgal_polyhedron_closest_d(ptree, xq,yq,zq, xn,yn,zn)
type(c_ptr),intent(in) :: ptree
real(c_double),intent(in) :: xq,yq,zq
real(c_double),intent(out) :: xn,yn,zn
type(d3) :: query, near
query = d3(xq,yq,zq)
call polyhedron_closest(ptree, query, near)
xn = near%x
yn = near%y
zn = near%z
end subroutine
subroutine cgal_polyhedron_bbox_s(ptree, xmin,ymin,zmin, xmax,ymax,zmax)
type(c_ptr),intent(in) :: ptree
real(c_float),intent(out) :: xmin,ymin,zmin
real(c_float),intent(out) :: xmax,ymax,zmax
type(d3) :: min, max
call polyhedron_bbox(ptree, min, max)
xmin = real(min%x,c_float)
ymin = real(min%y,c_float)
zmin = real(min%z,c_float)
xmax = real(max%x,c_float)
ymax = real(max%y,c_float)
zmax = real(max%z,c_float)
end subroutine
subroutine cgal_polyhedron_bbox_d(ptree, xmin,ymin,zmin, xmax,ymax,zmax)
type(c_ptr),intent(in) :: ptree
real(c_double),intent(out) :: xmin,ymin,zmin
real(c_double),intent(out) :: xmax,ymax,zmax
type(d3) :: min, max
call polyhedron_bbox(ptree, min, max)
xmin = min%x
ymin = min%y
zmin = min%z
xmax = max%x
ymax = max%y
zmax = max%z
end subroutine
function cgal_polyhedron_intersects_ray_s(ptree, xc,yc,zc, a,b,c) result(res)
logical :: res
type(c_ptr),intent(in) :: ptree
real(c_float),intent(in) :: xc,yc,zc, a,b,c
type(d3) :: origin, vec
origin = d3(real(xc,c_double),real(yc,c_double),real(zc,c_double))
vec = d3(real(a,c_double),real(b,c_double),real(c,c_double))
res = polyhedron_intersects_ray(ptree, origin, vec)
end function
function cgal_polyhedron_intersects_ray_d(ptree, xc,yc,zc, a,b,c) result(res)
logical :: res
type(c_ptr),intent(in) :: ptree
real(c_double),intent(in) :: xc,yc,zc, a,b,c
type(d3) :: origin, vec
origin = d3(xc,yc,zc)
vec = d3(a,b,c)
res = polyhedron_intersects_ray(ptree, origin, vec)
end function
subroutine cgal_polyhedron_finalize(ptree)
type(c_ptr),intent(inout) :: ptree
call polyhedron_finalize(ptree)
end subroutine
end module CGAL_Polyhedra