/
image.py
410 lines (351 loc) · 14.7 KB
/
image.py
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
# pyresample, Resampling of remote sensing image data in python
#
# Copyright (C) 2010, 2015 Esben S. Nielsen
#
# This program is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# This program is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU General Public License for more details.
#
# You should have received a copy of the GNU General Public License
# along with this program. If not, see <http://www.gnu.org/licenses/>.
"""Handles resampling of images with assigned geometry definitions."""
from __future__ import absolute_import
import warnings
import numpy as np
from pyresample import geometry, grid, kd_tree, bilinear
class ImageContainer(object):
"""Holds image with geometry definition. Allows indexing with linesample
arrays.
Parameters
----------
image_data : numpy array
Image data
geo_def : object
Geometry definition
fill_value : int or None, optional
Set undetermined pixels to this value.
If fill_value is None a masked array is returned
with undetermined pixels masked
nprocs : int, optional
Number of processor cores to be used
Attributes
----------
image_data : numpy array
Image data
geo_def : object
Geometry definition
fill_value : int or None
Resample result fill value
nprocs : int
Number of processor cores to be used for geometry operations
"""
def __init__(self, image_data, geo_def, fill_value=0, nprocs=1):
"""Initialize ImageContainer."""
warnings.warn(
"Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead",
FutureWarning)
if type(geo_def).__name__ == "DynamicAreaDefinition":
geo_def = geo_def.freeze()
if not isinstance(image_data, (np.ndarray, np.ma.core.MaskedArray)):
raise TypeError('image_data must be either an ndarray'
' or a masked array')
elif ((image_data.ndim > geo_def.ndim + 1) or
(image_data.ndim < geo_def.ndim)):
raise ValueError(('Unexpected number of dimensions for '
'image_data: %s') % image_data.ndim)
for i, size in enumerate(geo_def.shape):
if image_data.shape[i] != size:
raise ValueError(('Size mismatch for image_data. Expected '
'size %s for dimension %s and got %s') %
(size, i, image_data.shape[i]))
self.shape = geo_def.shape
self.size = geo_def.size
self.ndim = geo_def.ndim
self.image_data = image_data
if image_data.ndim > geo_def.ndim:
self.channels = image_data.shape[-1]
else:
self.channels = 1
self.geo_def = geo_def
self.fill_value = fill_value
self.nprocs = nprocs
def __str__(self):
return 'Image:\n %s' % self.image_data.__str__()
def __repr__(self):
return self.image_data.__repr__()
def resample(self, target_geo_def):
"""Base method for resampling."""
raise NotImplementedError('Method "resample" is not implemented '
'in class %s' % self.__class__.__name__)
def get_array_from_linesample(self, row_indices, col_indices):
"""Samples from image based on index arrays.
Parameters
----------
row_indices : numpy array
Row indices. Dimensions must match col_indices
col_indices : numpy array
Col indices. Dimensions must match row_indices
Returns
-------
image_data : numpy_array
Resampled image data
"""
if self.geo_def.ndim != 2:
raise TypeError('Resampling from linesamples only makes sense '
'on 2D data')
return grid.get_image_from_linesample(row_indices, col_indices,
self.image_data,
self.fill_value)
def get_array_from_neighbour_info(self, *args, **kwargs):
"""Base method for resampling from preprocessed data."""
raise NotImplementedError('Method "get_array_from_neighbour_info" is '
'not implemented in class %s' %
self.__class__.__name__)
class ImageContainerQuick(ImageContainer):
"""Holds image with area definition. ' Allows quick resampling within area.
Parameters
----------
image_data : numpy array
Image data
geo_def : object
Area definition as AreaDefinition object
fill_value : int or None, optional
Set undetermined pixels to this value.
If fill_value is None a masked array is returned
with undetermined pixels masked
nprocs : int, optional
Number of processor cores to be used for geometry operations
segments : int or None
Number of segments to use when resampling.
If set to None an estimate will be calculated
Attributes
----------
image_data : numpy array
Image data
geo_def : object
Area definition as AreaDefinition object
fill_value : int or None
Resample result fill value
If fill_value is None a masked array is returned
with undetermined pixels masked
nprocs : int
Number of processor cores to be used
segments : int or None
Number of segments to use when resampling
"""
def __init__(self, image_data, geo_def, fill_value=0, nprocs=1,
segments=None):
if not isinstance(geo_def, geometry.AreaDefinition):
raise TypeError('area_def must be of type '
'geometry.AreaDefinition')
super(ImageContainerQuick, self).__init__(image_data, geo_def,
fill_value=fill_value,
nprocs=nprocs)
self.segments = segments
def resample(self, target_area_def):
"""Resamples image to area definition using nearest neighbour approach
in projection coordinates.
Parameters
----------
target_area_def : object
Target area definition as AreaDefinition object
Returns
-------
image_container : object
ImageContainerQuick object of resampled area
"""
resampled_image = grid.get_resampled_image(target_area_def,
self.geo_def,
self.image_data,
fill_value=self.fill_value,
nprocs=self.nprocs,
segments=self.segments)
return ImageContainerQuick(resampled_image, target_area_def,
fill_value=self.fill_value,
nprocs=self.nprocs, segments=self.segments)
class ImageContainerNearest(ImageContainer):
"""Holds image with geometry definition. Allows nearest neighbour to new
geometry definition.
Parameters
----------
image_data : numpy array
Image data
geo_def : object
Geometry definition
radius_of_influence : float
Cut off distance in meters
epsilon : float, optional
Allowed uncertainty in meters. Increasing uncertainty
reduces execution time
fill_value : int or None, optional
Set undetermined pixels to this value.
If fill_value is None a masked array is returned
with undetermined pixels masked
reduce_data : bool, optional
Perform coarse data reduction before resampling in order
to reduce execution time
nprocs : int, optional
Number of processor cores to be used for geometry operations
segments : int or None
Number of segments to use when resampling.
If set to None an estimate will be calculated
Attributes
----------
image_data : numpy array
Image data
geo_def : object
Geometry definition
radius_of_influence : float
Cut off distance in meters
epsilon : float
Allowed uncertainty in meters
fill_value : int or None
Resample result fill value
reduce_data : bool
Perform coarse data reduction before resampling
nprocs : int
Number of processor cores to be used
segments : int or None
Number of segments to use when resampling
"""
def __init__(self, image_data, geo_def, radius_of_influence, epsilon=0,
fill_value=0, reduce_data=True, nprocs=1, segments=None):
super(ImageContainerNearest, self).__init__(image_data, geo_def,
fill_value=fill_value,
nprocs=nprocs)
self.radius_of_influence = radius_of_influence
self.epsilon = epsilon
self.reduce_data = reduce_data
self.segments = segments
def resample(self, target_geo_def):
"""Resamples image to area definition using nearest neighbour approach.
Parameters
----------
target_geo_def : object
Target geometry definition
Returns
-------
image_container : object
ImageContainerNearest object of resampled geometry
"""
if self.image_data.ndim > 2 and self.ndim > 1:
image_data = self.image_data.reshape(self.image_data.shape[0] *
self.image_data.shape[1],
self.image_data.shape[2])
else:
image_data = self.image_data.ravel()
resampled_image = \
kd_tree.resample_nearest(self.geo_def,
image_data,
target_geo_def,
self.radius_of_influence,
epsilon=self.epsilon,
fill_value=self.fill_value,
nprocs=self.nprocs,
reduce_data=self.reduce_data,
segments=self.segments)
return ImageContainerNearest(resampled_image, target_geo_def,
self.radius_of_influence,
epsilon=self.epsilon,
fill_value=self.fill_value,
reduce_data=self.reduce_data,
nprocs=self.nprocs,
segments=self.segments)
class ImageContainerBilinear(ImageContainer):
"""Holds image with geometry definition. Allows bilinear to new geometry
definition.
Parameters
----------
image_data : numpy array
Image data
geo_def : object
Geometry definition
radius_of_influence : float
Cut off distance in meters
epsilon : float, optional
Allowed uncertainty in meters. Increasing uncertainty
reduces execution time
fill_value : int or None, optional
Set undetermined pixels to this value.
If fill_value is None a masked array is returned
with undetermined pixels masked
reduce_data : bool, optional
Perform coarse data reduction before resampling in order
to reduce execution time
nprocs : int, optional
Number of processor cores to be used for geometry operations
segments : int or None
Number of segments to use when resampling.
If set to None an estimate will be calculated
Attributes
----------
image_data : numpy array
Image data
geo_def : object
Geometry definition
radius_of_influence : float
Cut off distance in meters
epsilon : float
Allowed uncertainty in meters
fill_value : int or None
Resample result fill value
reduce_data : bool
Perform coarse data reduction before resampling
nprocs : int
Number of processor cores to be used
segments : int or None
Number of segments to use when resampling
"""
def __init__(self, image_data, geo_def, radius_of_influence, epsilon=0,
fill_value=0, reduce_data=False, nprocs=1, segments=None,
neighbours=32):
super(ImageContainerBilinear, self).__init__(image_data, geo_def,
fill_value=fill_value,
nprocs=nprocs)
self.radius_of_influence = radius_of_influence
self.epsilon = epsilon
self.reduce_data = reduce_data
self.segments = segments
self.neighbours = neighbours
def resample(self, target_geo_def):
"""Resamples image to area definition using bilinear approach.
Parameters
----------
target_geo_def : object
Target geometry definition
Returns
-------
image_container : object
ImageContainerBilinear object of resampled geometry
"""
image_data = self.image_data
try:
mask = image_data.mask.copy()
image_data = image_data.data.copy()
image_data[mask] = np.nan
except AttributeError:
pass
resampled_image = \
bilinear.resample_bilinear(image_data,
self.geo_def,
target_geo_def,
radius=self.radius_of_influence,
neighbours=self.neighbours,
epsilon=self.epsilon,
fill_value=self.fill_value,
nprocs=self.nprocs,
reduce_data=self.reduce_data,
segments=self.segments)
return ImageContainerBilinear(resampled_image, target_geo_def,
self.radius_of_influence,
epsilon=self.epsilon,
fill_value=self.fill_value,
reduce_data=self.reduce_data,
nprocs=self.nprocs,
segments=self.segments)