forked from GeoStat-Framework/GSTools
-
Notifications
You must be signed in to change notification settings - Fork 0
/
estimator.pyx
314 lines (272 loc) · 9.4 KB
/
estimator.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
#cython: language_level=3, boundscheck=False, wraparound=False, cdivision=True
# distutils: language = c++
# -*- coding: utf-8 -*-
"""
This is the variogram estimater, implemented in cython.
"""
import numpy as np
cimport cython
from cython.parallel import prange, parallel
from libcpp.vector cimport vector
from libc.math cimport fabs, sqrt, atan2, acos, asin, sin, cos, M_PI
cimport numpy as np
DTYPE = np.double
ctypedef np.double_t DTYPE_t
cdef inline double _distance_1d(
const double[:] x,
const double[:] y,
const double[:] z,
const int i,
const int j
) nogil:
return sqrt((x[i] - x[j]) * (x[i] - x[j]))
cdef inline double _distance_2d(
const double[:] x,
const double[:] y,
const double[:] z,
const int i,
const int j
) nogil:
return sqrt((x[i] - x[j]) * (x[i] - x[j]) + (y[i] - y[j]) * (y[i] - y[j]))
cdef inline double _distance_3d(
const double[:] x,
const double[:] y,
const double[:] z,
const int i,
const int j
) nogil:
return sqrt((x[i] - x[j]) * (x[i] - x[j]) +
(y[i] - y[j]) * (y[i] - y[j]) +
(z[i] - z[j]) * (z[i] - z[j]))
cdef inline bint _angle_test_1d(
const double[:] x,
const double[:] y,
const double[:] z,
const double[:] angles,
const double angles_tol,
const int i,
const int j
) nogil:
return True
cdef inline bint _angle_test_2d(
const double[:] x,
const double[:] y,
const double[:] z,
const double[:] angles,
const double angles_tol,
const int i,
const int j
) nogil:
cdef double dx = x[i] - x[j]
cdef double dy = y[i] - y[j]
# azimuth
cdef double phi1 = atan2(dy,dx) % (2.0 * M_PI)
cdef double phi2 = atan2(-dy,-dx) % (2.0 * M_PI)
# check both directions (+/-)
cdef bint dir1 = fabs(phi1 - angles[0]) <= angles_tol
cdef bint dir2 = fabs(phi2 - angles[0]) <= angles_tol
return dir1 or dir2
cdef inline bint _angle_test_3d(
const double[:] x,
const double[:] y,
const double[:] z,
const double[:] angles,
const double angles_tol,
const int i,
const int j
) nogil:
cdef double dx = x[i] - x[j]
cdef double dy = y[i] - y[j]
cdef double dz = z[i] - z[j]
cdef double dr = sqrt(dx**2 + dy**2 + dz**2)
# azimuth
cdef double phi1 = atan2(dy, dx) % (2.0 * M_PI)
cdef double phi2 = atan2(-dy, -dx) % (2.0 * M_PI)
# elevation
cdef double theta1 = acos(dz / dr)
cdef double theta2 = acos(-dz / dr)
# definitions for great-circle distance (for tolerance check)
cdef double dx1 = sin(theta1) * cos(phi1) - sin(angles[1]) * cos(angles[0])
cdef double dy1 = sin(theta1) * sin(phi1) - sin(angles[1]) * sin(angles[0])
cdef double dz1 = cos(theta1) - cos(angles[1])
cdef double dx2 = sin(theta2) * cos(phi2) - sin(angles[1]) * cos(angles[0])
cdef double dy2 = sin(theta2) * sin(phi2) - sin(angles[1]) * sin(angles[0])
cdef double dz2 = cos(theta2) - cos(angles[1])
cdef double dist1 = 2.0 * asin(sqrt(dx1**2 + dy1**2 + dz1**2) * 0.5)
cdef double dist2 = 2.0 * asin(sqrt(dx2**2 + dy2**2 + dz2**2) * 0.5)
# check both directions (+/-)
cdef bint dir1 = dist1 <= angles_tol
cdef bint dir2 = dist2 <= angles_tol
return dir1 or dir2
cdef inline double estimator_matheron(const double f_diff) nogil:
return f_diff * f_diff
cdef inline double estimator_cressie(const double f_diff) nogil:
return sqrt(fabs(f_diff))
ctypedef double (*_estimator_func)(const double) nogil
cdef inline void normalization_matheron(
vector[double]& variogram,
vector[long]& counts
):
cdef int i
for i in range(variogram.size()):
# avoid division by zero
if counts[i] == 0:
counts[i] = 1
variogram[i] /= (2. * counts[i])
cdef inline void normalization_cressie(
vector[double]& variogram,
vector[long]& counts
):
cdef int i
for i in range(variogram.size()):
# avoid division by zero
if counts[i] == 0:
counts[i] = 1
variogram[i] = (
(1./counts[i] * variogram[i])**4 /
(0.457 + 0.494 / counts[i] + 0.045 / counts[i]**2)
)
ctypedef void (*_normalization_func)(
vector[double]&,
vector[long]&
)
cdef _estimator_func choose_estimator_func(str estimator_type):
cdef _estimator_func estimator_func
if estimator_type == 'm':
estimator_func = estimator_matheron
elif estimator_type == 'c':
estimator_func = estimator_cressie
return estimator_func
cdef _normalization_func choose_estimator_normalization(str estimator_type):
cdef _normalization_func normalization_func
if estimator_type == 'm':
normalization_func = normalization_matheron
elif estimator_type == 'c':
normalization_func = normalization_cressie
return normalization_func
ctypedef double (*_dist_func)(
const double[:],
const double[:],
const double[:],
const int,
const int
) nogil
ctypedef bint (*_angle_test_func)(
const double[:],
const double[:],
const double[:],
const double[:],
const double,
const int,
const int
) nogil
def unstructured(
const double[:] f,
const double[:] bin_edges,
const double[:] x,
const double[:] y=None,
const double[:] z=None,
const double[:] angles=None,
const double angles_tol=0.436332,
str estimator_type='m'
):
if x.shape[0] != f.shape[0]:
raise ValueError('len(x) = {0} != len(f) = {1} '.
format(x.shape[0], f.shape[0]))
if bin_edges.shape[0] < 2:
raise ValueError('len(bin_edges) too small')
cdef _dist_func distance
cdef _angle_test_func angle_test
# 3d
if z is not None:
if z.shape[0] != f.shape[0]:
raise ValueError('len(z) = {0} != len(f) = {1} '.
format(z.shape[0], f.shape[0]))
distance = _distance_3d
angle_test = _angle_test_3d
# 2d
elif y is not None:
if y.shape[0] != f.shape[0]:
raise ValueError('len(y) = {0} != len(f) = {1} '.
format(y.shape[0], f.shape[0]))
distance = _distance_2d
angle_test = _angle_test_2d
# 1d
else:
distance = _distance_1d
angle_test = _angle_test_1d
if angles is not None:
if z is not None and angles.size < 2:
raise ValueError('3d requested but only one angle given')
if y is not None and angles.size < 1:
raise ValueError('2d with angle requested but no angle given')
if angles_tol <= 0:
raise ValueError('tolerance for angle search masks must be > 0')
cdef _estimator_func estimator_func = choose_estimator_func(estimator_type)
cdef _normalization_func normalization_func = (
choose_estimator_normalization(estimator_type)
)
cdef int i_max = bin_edges.shape[0] - 1
cdef int j_max = x.shape[0] - 1
cdef int k_max = x.shape[0]
cdef vector[double] variogram = vector[double](len(bin_edges)-1, 0.0)
cdef vector[long] counts = vector[long](len(bin_edges)-1, 0)
cdef int i, j, k
cdef DTYPE_t dist
for i in prange(i_max, nogil=True):
for j in range(j_max):
for k in range(j+1, k_max):
dist = distance(x, y, z, k, j)
if dist >= bin_edges[i] and dist < bin_edges[i+1]:
if angles is None or angle_test(x, y, z, angles, angles_tol, k, j):
counts[i] += 1
variogram[i] += estimator_func(f[k] - f[j])
normalization_func(variogram, counts)
return np.asarray(variogram), np.asarray(counts)
def structured(const double[:,:,:] f, str estimator_type='m'):
cdef _estimator_func estimator_func = choose_estimator_func(estimator_type)
cdef _normalization_func normalization_func = (
choose_estimator_normalization(estimator_type)
)
cdef int i_max = f.shape[0] - 1
cdef int j_max = f.shape[1]
cdef int k_max = f.shape[2]
cdef int l_max = i_max + 1
cdef vector[double] variogram = vector[double](l_max, 0.0)
cdef vector[long] counts = vector[long](l_max, 0)
cdef int i, j, k, l
with nogil, parallel():
for i in range(i_max):
for j in range(j_max):
for k in range(k_max):
for l in prange(1, l_max-i):
counts[l] += 1
variogram[l] += estimator_func(f[i,j,k] - f[i+l,j,k])
normalization_func(variogram, counts)
return np.asarray(variogram)
def ma_structured(
const double[:,:,:] f,
const bint[:,:,:] mask,
str estimator_type='m'
):
cdef _estimator_func estimator_func = choose_estimator_func(estimator_type)
cdef _normalization_func normalization_func = (
choose_estimator_normalization(estimator_type)
)
cdef int i_max = f.shape[0] - 1
cdef int j_max = f.shape[1]
cdef int k_max = f.shape[2]
cdef int l_max = i_max + 1
cdef vector[double] variogram = vector[double](l_max, 0.0)
cdef vector[long] counts = vector[long](l_max, 0)
cdef int i, j, k, l
with nogil, parallel():
for i in range(i_max):
for j in range(j_max):
for k in range(k_max):
for l in prange(1, l_max-i):
if not mask[i,j,k] and not mask[i+l,j,k]:
counts[l] += 1
variogram[l] += estimator_func(f[i,j,k] - f[i+l,j,k])
normalization_func(variogram, counts)
return np.asarray(variogram)