/
resample.py
1467 lines (1209 loc) · 56.9 KB
/
resample.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
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
#!/usr/bin/env python
# -*- coding: utf-8 -*-
# Copyright (c) 2015-2018 Satpy developers
#
# This file is part of satpy.
#
# satpy 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.
#
# satpy 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
# satpy. If not, see <http://www.gnu.org/licenses/>.
"""Resampling in Satpy.
Satpy provides multiple resampling algorithms for resampling geolocated
data to uniform projected grids. The easiest way to perform resampling in
Satpy is through the :class:`~satpy.scene.Scene` object's
:meth:`~satpy.scene.Scene.resample` method. Additional utility functions are
also available to assist in resampling data. Below is more information on
resampling with Satpy as well as links to the relevant API documentation for
available keyword arguments.
Resampling algorithms
---------------------
.. csv-table:: Available Resampling Algorithms
:header-rows: 1
:align: center
"Resampler", "Description", "Related"
"nearest", "Nearest Neighbor", :class:`~satpy.resample.KDTreeResampler`
"ewa", "Elliptical Weighted Averaging", :class:`~pyresample.ewa.DaskEWAResampler`
"ewa_legacy", "Elliptical Weighted Averaging (Legacy)", :class:`~pyresample.ewa.LegacyDaskEWAResampler`
"native", "Native", :class:`~satpy.resample.NativeResampler`
"bilinear", "Bilinear", :class:`~satpy.resample.BilinearResampler`
"bucket_avg", "Average Bucket Resampling", :class:`~satpy.resample.BucketAvg`
"bucket_sum", "Sum Bucket Resampling", :class:`~satpy.resample.BucketSum`
"bucket_count", "Count Bucket Resampling", :class:`~satpy.resample.BucketCount`
"bucket_fraction", "Fraction Bucket Resampling", :class:`~satpy.resample.BucketFraction`
"gradient_search", "Gradient Search Resampling", :class:`~pyresample.gradient.GradientSearchResampler`
The resampling algorithm used can be specified with the ``resampler`` keyword
argument and defaults to ``nearest``:
.. code-block:: python
>>> scn = Scene(...)
>>> euro_scn = scn.resample('euro4', resampler='nearest')
.. warning::
Some resampling algorithms expect certain forms of data. For example, the
EWA resampling expects polar-orbiting swath data and prefers if the data
can be broken in to "scan lines". See the API documentation for a specific
algorithm for more information.
Resampling for comparison and composites
----------------------------------------
While all the resamplers can be used to put datasets of different resolutions
on to a common area, the 'native' resampler is designed to match datasets to
one resolution in the dataset's original projection. This is extremely useful
when generating composites between bands of different resolutions.
.. code-block:: python
>>> new_scn = scn.resample(resampler='native')
By default this resamples to the
:meth:`highest resolution area <satpy.scene.Scene.finest_area>` (smallest footprint per
pixel) shared between the loaded datasets. You can easily specify the lowest
resolution area:
.. code-block:: python
>>> new_scn = scn.resample(scn.coarsest_area(), resampler='native')
Providing an area that is neither the minimum or maximum resolution area
may work, but behavior is currently undefined.
Caching for geostationary data
------------------------------
Satpy will do its best to reuse calculations performed to resample datasets,
but it can only do this for the current processing and will lose this
information when the process/script ends. Some resampling algorithms, like
``nearest`` and ``bilinear``, can benefit by caching intermediate data on disk in the directory
specified by `cache_dir` and using it next time. This is most beneficial with
geostationary satellite data where the locations of the source data and the
target pixels don't change over time.
>>> new_scn = scn.resample('euro4', cache_dir='/path/to/cache_dir')
See the documentation for specific algorithms to see availability and
limitations of caching for that algorithm.
Create custom area definition
-----------------------------
See :class:`pyresample.geometry.AreaDefinition` for information on creating
areas that can be passed to the resample method::
>>> from pyresample.geometry import AreaDefinition
>>> my_area = AreaDefinition(...)
>>> local_scene = scn.resample(my_area)
Create dynamic area definition
------------------------------
See :class:`pyresample.geometry.DynamicAreaDefinition` for more information.
Examples coming soon...
Store area definitions
----------------------
Area definitions can be saved to a custom YAML file (see
`pyresample's writing to disk <http://pyresample.readthedocs.io/en/stable/geometry_utils.html#writing-to-disk>`_)
and loaded using pyresample's utility methods
(`pyresample's loading from disk <http://pyresample.readthedocs.io/en/stable/geometry_utils.html#loading-from-disk>`_)::
>>> from pyresample import load_area
>>> my_area = load_area('my_areas.yaml', 'my_area')
Or using :func:`satpy.resample.get_area_def`, which will search through all
``areas.yaml`` files in your ``SATPY_CONFIG_PATH``::
>>> from satpy.resample import get_area_def
>>> area_eurol = get_area_def("eurol")
For examples of area definitions, see the file ``etc/areas.yaml`` that is
included with Satpy and where all the area definitions shipped with Satpy are
defined.
"""
import hashlib
import json
import os
import warnings
from logging import getLogger
from weakref import WeakValueDictionary
import dask
import dask.array as da
import numpy as np
import pyresample
import xarray as xr
import zarr
from packaging import version
from pyresample.ewa import fornav, ll2cr
from pyresample.geometry import SwathDefinition
from satpy.utils import PerformanceWarning, get_legacy_chunk_size
try:
from math import lcm # type: ignore
except ImportError:
def lcm(a, b):
"""Get 'Least Common Multiple' with Python 3.8 compatibility."""
from math import gcd
return abs(a * b) // gcd(a, b)
try:
from pyresample.resampler import BaseResampler as PRBaseResampler
except ImportError:
PRBaseResampler = None
try:
from pyresample.gradient import GradientSearchResampler
except ImportError:
GradientSearchResampler = None
try:
from pyresample.ewa import DaskEWAResampler, LegacyDaskEWAResampler
except ImportError:
DaskEWAResampler = LegacyDaskEWAResampler = None
from satpy._config import config_search_paths, get_config_path
LOG = getLogger(__name__)
CHUNK_SIZE = get_legacy_chunk_size()
CACHE_SIZE = 10
NN_COORDINATES = {'valid_input_index': ('y1', 'x1'),
'valid_output_index': ('y2', 'x2'),
'index_array': ('y2', 'x2', 'z2')}
BIL_COORDINATES = {'bilinear_s': ('x1', ),
'bilinear_t': ('x1', ),
'slices_x': ('x1', 'n'),
'slices_y': ('x1', 'n'),
'mask_slices': ('x1', 'n'),
'out_coords_x': ('x2', ),
'out_coords_y': ('y2', )}
resamplers_cache: "WeakValueDictionary[tuple, object]" = WeakValueDictionary()
PR_USE_SKIPNA = version.parse(pyresample.__version__) > version.parse("1.17.0")
def hash_dict(the_dict, the_hash=None):
"""Calculate a hash for a dictionary."""
if the_hash is None:
the_hash = hashlib.sha1() # nosec
the_hash.update(json.dumps(the_dict, sort_keys=True).encode('utf-8'))
return the_hash
def get_area_file():
"""Find area file(s) to use.
The files are to be named `areas.yaml` or `areas.def`.
"""
paths = config_search_paths('areas.yaml')
if paths:
return paths
else:
return get_config_path('areas.def')
def get_area_def(area_name):
"""Get the definition of *area_name* from file.
The file is defined to use is to be placed in the $SATPY_CONFIG_PATH
directory, and its name is defined in satpy's configuration file.
"""
try:
from pyresample import parse_area_file
except ImportError:
from pyresample.utils import parse_area_file
return parse_area_file(get_area_file(), area_name)[0]
def add_xy_coords(data_arr, area, crs=None):
"""Assign x/y coordinates to DataArray from provided area.
If 'x' and 'y' coordinates already exist then they will not be added.
Args:
data_arr (xarray.DataArray): data object to add x/y coordinates to
area (pyresample.geometry.AreaDefinition): area providing the
coordinate data.
crs (pyproj.crs.CRS or None): CRS providing additional information
about the area's coordinate reference system if available.
Requires pyproj 2.0+.
Returns (xarray.DataArray): Updated DataArray object
"""
if 'x' in data_arr.coords and 'y' in data_arr.coords:
# x/y coords already provided
return data_arr
if 'x' not in data_arr.dims or 'y' not in data_arr.dims:
# no defined x and y dimensions
return data_arr
if not hasattr(area, 'get_proj_vectors'):
return data_arr
x, y = area.get_proj_vectors()
# convert to DataArrays
y_attrs = {}
x_attrs = {}
if crs is not None:
units = crs.axis_info[0].unit_name
# fix udunits/CF standard units
units = units.replace('metre', 'meter')
if units == 'degree':
y_attrs['units'] = 'degrees_north'
x_attrs['units'] = 'degrees_east'
else:
y_attrs['units'] = units
x_attrs['units'] = units
y = xr.DataArray(y, dims=('y',), attrs=y_attrs)
x = xr.DataArray(x, dims=('x',), attrs=x_attrs)
return data_arr.assign_coords(y=y, x=x)
def add_crs_xy_coords(data_arr, area):
"""Add :class:`pyproj.crs.CRS` and x/y or lons/lats to coordinates.
For SwathDefinition or GridDefinition areas this will add a
`crs` coordinate and coordinates for the 2D arrays of `lons` and `lats`.
For AreaDefinition areas this will add a `crs` coordinate and the
1-dimensional `x` and `y` coordinate variables.
Args:
data_arr (xarray.DataArray): DataArray to add the 'crs'
coordinate.
area (pyresample.geometry.AreaDefinition): Area to get CRS
information from.
"""
# add CRS object if pyproj 2.0+
try:
from pyproj import CRS
except ImportError:
LOG.debug("Could not add 'crs' coordinate with pyproj<2.0")
crs = None
else:
# default lat/lon projection
latlon_proj = "+proj=latlong +datum=WGS84 +ellps=WGS84"
# otherwise get it from the area definition
if hasattr(area, 'crs'):
crs = area.crs
else:
proj_str = getattr(area, 'proj_str', latlon_proj)
crs = CRS.from_string(proj_str)
data_arr = data_arr.assign_coords(crs=crs)
# Add x/y coordinates if possible
if isinstance(area, SwathDefinition):
# add lon/lat arrays for swath definitions
# SwathDefinitions created by Satpy should be assigning DataArray
# objects as the lons/lats attributes so use those directly to
# maintain original .attrs metadata (instead of converting to dask
# array).
lons = area.lons
lats = area.lats
lons.attrs.setdefault('standard_name', 'longitude')
lons.attrs.setdefault('long_name', 'longitude')
lons.attrs.setdefault('units', 'degrees_east')
lats.attrs.setdefault('standard_name', 'latitude')
lats.attrs.setdefault('long_name', 'latitude')
lats.attrs.setdefault('units', 'degrees_north')
# See https://github.com/pydata/xarray/issues/3068
# data_arr = data_arr.assign_coords(longitude=lons, latitude=lats)
else:
# Gridded data (AreaDefinition/StackedAreaDefinition)
data_arr = add_xy_coords(data_arr, area, crs=crs)
return data_arr
def update_resampled_coords(old_data, new_data, new_area):
"""Add coordinate information to newly resampled DataArray.
Args:
old_data (xarray.DataArray): Old data before resampling.
new_data (xarray.DataArray): New data after resampling.
new_area (pyresample.geometry.BaseDefinition): Area definition
for the newly resampled data.
"""
# copy over other non-x/y coordinates
# this *MUST* happen before we set 'crs' below otherwise any 'crs'
# coordinate in the coordinate variables we are copying will overwrite the
# 'crs' coordinate we just assigned to the data
ignore_coords = ('y', 'x', 'crs')
new_coords = {}
for cname, cval in old_data.coords.items():
# we don't want coordinates that depended on the old x/y dimensions
has_ignored_dims = any(dim in cval.dims for dim in ignore_coords)
if cname in ignore_coords or has_ignored_dims:
continue
new_coords[cname] = cval
new_data = new_data.assign_coords(**new_coords)
# add crs, x, and y coordinates
new_data = add_crs_xy_coords(new_data, new_area)
return new_data
class BaseResampler(object):
"""Base abstract resampler class."""
def __init__(self, source_geo_def, target_geo_def):
"""Initialize resampler with geolocation information.
Args:
source_geo_def (SwathDefinition, AreaDefinition):
Geolocation definition for the data to be resampled
target_geo_def (CoordinateDefinition, AreaDefinition):
Geolocation definition for the area to resample data to.
"""
self.source_geo_def = source_geo_def
self.target_geo_def = target_geo_def
def get_hash(self, source_geo_def=None, target_geo_def=None, **kwargs):
"""Get hash for the current resample with the given *kwargs*."""
if source_geo_def is None:
source_geo_def = self.source_geo_def
if target_geo_def is None:
target_geo_def = self.target_geo_def
the_hash = source_geo_def.update_hash()
target_geo_def.update_hash(the_hash)
hash_dict(kwargs, the_hash)
return the_hash.hexdigest()
def precompute(self, **kwargs):
"""Do the precomputation.
This is an optional step if the subclass wants to implement more
complex features like caching or can share some calculations
between multiple datasets to be processed.
"""
return None
def compute(self, data, **kwargs):
"""Do the actual resampling.
This must be implemented by subclasses.
"""
raise NotImplementedError
def resample(self, data, cache_dir=None, mask_area=None, **kwargs):
"""Resample `data` by calling `precompute` and `compute` methods.
Only certain resampling classes may use `cache_dir` and the `mask`
provided when `mask_area` is True. The return value of calling the
`precompute` method is passed as the `cache_id` keyword argument
of the `compute` method, but may not be used directly for caching. It
is up to the individual resampler subclasses to determine how this
is used.
Args:
data (xarray.DataArray): Data to be resampled
cache_dir (str): directory to cache precomputed results
(default False, optional)
mask_area (bool): Mask geolocation data where data values are
invalid. This should be used when data values
may affect what neighbors are considered valid.
Returns (xarray.DataArray): Data resampled to the target area
"""
# default is to mask areas for SwathDefinitions
if mask_area is None and isinstance(
self.source_geo_def, SwathDefinition):
mask_area = True
if mask_area:
if isinstance(self.source_geo_def, SwathDefinition):
geo_dims = self.source_geo_def.lons.dims
else:
geo_dims = ('y', 'x')
flat_dims = [dim for dim in data.dims if dim not in geo_dims]
if np.issubdtype(data.dtype, np.integer):
kwargs['mask'] = data == data.attrs.get('_FillValue', np.iinfo(data.dtype.type).max)
else:
kwargs['mask'] = data.isnull()
kwargs['mask'] = kwargs['mask'].all(dim=flat_dims)
cache_id = self.precompute(cache_dir=cache_dir, **kwargs)
return self.compute(data, cache_id=cache_id, **kwargs)
def _create_cache_filename(self, cache_dir, prefix='',
fmt='.zarr', **kwargs):
"""Create filename for the cached resampling parameters."""
hash_str = self.get_hash(**kwargs)
return os.path.join(cache_dir, prefix + hash_str + fmt)
class KDTreeResampler(BaseResampler):
"""Resample using a KDTree-based nearest neighbor algorithm.
This resampler implements on-disk caching when the `cache_dir` argument
is provided to the `resample` method. This should provide significant
performance improvements on consecutive resampling of geostationary data.
It is not recommended to provide `cache_dir` when the `mask` keyword
argument is provided to `precompute` which occurs by default for
`SwathDefinition` source areas.
Args:
cache_dir (str): Long term storage directory for intermediate
results.
mask (bool): Force resampled data's invalid pixel mask to be used
when searching for nearest neighbor pixels. By
default this is True for SwathDefinition source
areas and False for all other area definition types.
radius_of_influence (float): Search radius cut off distance in meters
epsilon (float): Allowed uncertainty in meters. Increasing uncertainty
reduces execution time.
"""
def __init__(self, source_geo_def, target_geo_def):
"""Init KDTreeResampler."""
super(KDTreeResampler, self).__init__(source_geo_def, target_geo_def)
self.resampler = None
self._index_caches = {}
def precompute(self, mask=None, radius_of_influence=None, epsilon=0,
cache_dir=None, **kwargs):
"""Create a KDTree structure and store it for later use.
Note: The `mask` keyword should be provided if geolocation may be valid
where data points are invalid.
"""
from pyresample.kd_tree import XArrayResamplerNN
del kwargs
if mask is not None and cache_dir is not None:
LOG.warning("Mask and cache_dir both provided to nearest "
"resampler. Cached parameters are affected by "
"masked pixels. Will not cache results.")
cache_dir = None
if radius_of_influence is None and not hasattr(self.source_geo_def, 'geocentric_resolution'):
radius_of_influence = self._adjust_radius_of_influence(radius_of_influence)
kwargs = dict(source_geo_def=self.source_geo_def,
target_geo_def=self.target_geo_def,
radius_of_influence=radius_of_influence,
neighbours=1,
epsilon=epsilon)
if self.resampler is None:
# FIXME: We need to move all of this caching logic to pyresample
self.resampler = XArrayResamplerNN(**kwargs)
try:
self.load_neighbour_info(cache_dir, mask=mask, **kwargs)
LOG.debug("Read pre-computed kd-tree parameters")
except IOError:
LOG.debug("Computing kd-tree parameters")
self.resampler.get_neighbour_info(mask=mask)
self.save_neighbour_info(cache_dir, mask=mask, **kwargs)
def _adjust_radius_of_influence(self, radius_of_influence):
"""Adjust radius of influence."""
warnings.warn(
"Upgrade 'pyresample' for a more accurate default 'radius_of_influence'.",
stacklevel=3
)
try:
radius_of_influence = self.source_geo_def.lons.resolution * 3
except AttributeError:
try:
radius_of_influence = max(abs(self.source_geo_def.pixel_size_x),
abs(self.source_geo_def.pixel_size_y)) * 3
except AttributeError:
radius_of_influence = 1000
except TypeError:
radius_of_influence = 10000
return radius_of_influence
def _apply_cached_index(self, val, idx_name, persist=False):
"""Reassign resampler index attributes."""
if isinstance(val, np.ndarray):
val = da.from_array(val, chunks=CHUNK_SIZE)
elif persist and isinstance(val, da.Array):
val = val.persist()
setattr(self.resampler, idx_name, val)
return val
def _check_numpy_cache(self, cache_dir, mask=None,
**kwargs):
"""Check if there's Numpy cache file and convert it to zarr."""
if cache_dir is None:
return
fname_np = self._create_cache_filename(cache_dir,
prefix='resample_lut-',
mask=mask, fmt='.npz',
**kwargs)
fname_zarr = self._create_cache_filename(cache_dir, prefix='nn_lut-',
mask=mask, fmt='.zarr',
**kwargs)
LOG.debug("Check if %s exists", fname_np)
if os.path.exists(fname_np) and not os.path.exists(fname_zarr):
import warnings
warnings.warn(
"Using Numpy files as resampling cache is deprecated.",
stacklevel=3
)
LOG.warning("Converting resampling LUT from .npz to .zarr")
zarr_out = xr.Dataset()
with np.load(fname_np, 'r') as fid:
for idx_name, coord in NN_COORDINATES.items():
zarr_out[idx_name] = (coord, fid[idx_name])
# Write indices to Zarr file
zarr_out.to_zarr(fname_zarr)
LOG.debug("Resampling LUT saved to %s", fname_zarr)
def load_neighbour_info(self, cache_dir, mask=None, **kwargs):
"""Read index arrays from either the in-memory or disk cache."""
mask_name = getattr(mask, 'name', None)
cached = {}
self._check_numpy_cache(cache_dir, mask=mask_name, **kwargs)
for idx_name in NN_COORDINATES:
if mask_name in self._index_caches:
cached[idx_name] = self._apply_cached_index(
self._index_caches[mask_name][idx_name], idx_name)
elif cache_dir:
try:
filename = self._create_cache_filename(
cache_dir, prefix='nn_lut-',
mask=mask_name, **kwargs)
fid = zarr.open(filename, 'r')
cache = np.array(fid[idx_name])
if idx_name == 'valid_input_index':
# valid input index array needs to be boolean
cache = cache.astype(bool)
except ValueError:
raise IOError
cache = self._apply_cached_index(cache, idx_name)
cached[idx_name] = cache
else:
raise IOError
self._index_caches[mask_name] = cached
def save_neighbour_info(self, cache_dir, mask=None, **kwargs):
"""Cache resampler's index arrays if there is a cache dir."""
if cache_dir:
mask_name = getattr(mask, 'name', None)
cache = self._read_resampler_attrs()
filename = self._create_cache_filename(
cache_dir, prefix='nn_lut-', mask=mask_name, **kwargs)
LOG.info('Saving kd_tree neighbour info to %s', filename)
zarr_out = xr.Dataset()
for idx_name, coord in NN_COORDINATES.items():
# update the cache in place with persisted dask arrays
cache[idx_name] = self._apply_cached_index(cache[idx_name],
idx_name,
persist=True)
zarr_out[idx_name] = (coord, cache[idx_name])
# Write indices to Zarr file
zarr_out.to_zarr(filename)
self._index_caches[mask_name] = cache
# Delete the kdtree, it's not needed anymore
self.resampler.delayed_kdtree = None
def _read_resampler_attrs(self):
"""Read certain attributes from the resampler for caching."""
return {attr_name: getattr(self.resampler, attr_name)
for attr_name in NN_COORDINATES}
def compute(self, data, weight_funcs=None, fill_value=np.nan,
with_uncert=False, **kwargs):
"""Resample data."""
del kwargs
LOG.debug("Resampling %s", str(data.name))
res = self.resampler.get_sample_from_neighbour_info(data, fill_value)
return update_resampled_coords(data, res, self.target_geo_def)
class _LegacySatpyEWAResampler(BaseResampler):
"""Resample using an elliptical weighted averaging algorithm.
This algorithm does **not** use caching or any externally provided data
mask (unlike the 'nearest' resampler).
This algorithm works under the assumption that the data is observed
one scan line at a time. However, good results can still be achieved
for non-scan based data provided `rows_per_scan` is set to the
number of rows in the entire swath or by setting it to `None`.
Args:
rows_per_scan (int, None):
Number of data rows for every observed scanline. If None then the
entire swath is treated as one large scanline.
weight_count (int):
number of elements to create in the gaussian weight table.
Default is 10000. Must be at least 2
weight_min (float):
the minimum value to store in the last position of the
weight table. Default is 0.01, which, with a
`weight_distance_max` of 1.0 produces a weight of 0.01
at a grid cell distance of 1.0. Must be greater than 0.
weight_distance_max (float):
distance in grid cell units at which to
apply a weight of `weight_min`. Default is
1.0. Must be greater than 0.
weight_delta_max (float):
maximum distance in grid cells in each grid
dimension over which to distribute a single swath cell.
Default is 10.0.
weight_sum_min (float):
minimum weight sum value. Cells whose weight sums
are less than `weight_sum_min` are set to the grid fill value.
Default is EPSILON.
maximum_weight_mode (bool):
If False (default), a weighted average of
all swath cells that map to a particular grid cell is used.
If True, the swath cell having the maximum weight of all
swath cells that map to a particular grid cell is used. This
option should be used for coded/category data, i.e. snow cover.
"""
def __init__(self, source_geo_def, target_geo_def):
"""Init _LegacySatpyEWAResampler."""
warnings.warn(
"A new version of pyresample is available. Please "
"upgrade to get access to a newer 'ewa' and "
"'ewa_legacy' resampler.",
stacklevel=2
)
super(_LegacySatpyEWAResampler, self).__init__(source_geo_def, target_geo_def)
self.cache = {}
def resample(self, *args, **kwargs):
"""Run precompute and compute methods.
.. note::
This sets the default of 'mask_area' to False since it is
not needed in EWA resampling currently.
"""
kwargs.setdefault('mask_area', False)
return super(_LegacySatpyEWAResampler, self).resample(*args, **kwargs)
def _call_ll2cr(self, lons, lats, target_geo_def, swath_usage=0):
"""Wrap ll2cr() for handling dask delayed calls better."""
new_src = SwathDefinition(lons, lats)
swath_points_in_grid, cols, rows = ll2cr(new_src, target_geo_def)
# FIXME: How do we check swath usage/coverage if we only do this
# per-block
# # Determine if enough of the input swath was used
# grid_name = getattr(self.target_geo_def, "name", "N/A")
# fraction_in = swath_points_in_grid / float(lons.size)
# swath_used = fraction_in > swath_usage
# if not swath_used:
# LOG.info("Data does not fit in grid %s because it only %f%% of "
# "the swath is used" %
# (grid_name, fraction_in * 100))
# raise RuntimeError("Data does not fit in grid %s" % (grid_name,))
# else:
# LOG.debug("Data fits in grid %s and uses %f%% of the swath",
# grid_name, fraction_in * 100)
return np.stack([cols, rows], axis=0)
def precompute(self, cache_dir=None, swath_usage=0, **kwargs):
"""Generate row and column arrays and store it for later use."""
if self.cache:
# this resampler should be used for one SwathDefinition
# no need to recompute ll2cr output again
return None
if kwargs.get('mask') is not None:
LOG.warning("'mask' parameter has no affect during EWA "
"resampling")
del kwargs
source_geo_def = self.source_geo_def
target_geo_def = self.target_geo_def
if cache_dir:
LOG.warning("'cache_dir' is not used by EWA resampling")
# Satpy/PyResample don't support dynamic grids out of the box yet
lons, lats = source_geo_def.get_lonlats()
if isinstance(lons, xr.DataArray):
# get dask arrays
lons = lons.data
lats = lats.data
# we are remapping to a static unchanging grid/area with all of
# its parameters specified
chunks = (2,) + lons.chunks
res = da.map_blocks(self._call_ll2cr, lons, lats,
target_geo_def, swath_usage,
dtype=lons.dtype, chunks=chunks, new_axis=[0])
cols = res[0]
rows = res[1]
# save the dask arrays in the class instance cache
# the on-disk cache will store the numpy arrays
self.cache = {
"rows": rows,
"cols": cols,
}
return None
def _call_fornav(self, cols, rows, target_geo_def, data,
grid_coverage=0, **kwargs):
"""Wrap fornav() to run as a dask delayed."""
num_valid_points, res = fornav(cols, rows, target_geo_def,
data, **kwargs)
if isinstance(data, tuple):
# convert 'res' from tuple of arrays to one array
res = np.stack(res)
num_valid_points = sum(num_valid_points)
grid_covered_ratio = num_valid_points / float(res.size)
grid_covered = grid_covered_ratio > grid_coverage
if not grid_covered:
msg = "EWA resampling only found %f%% of the grid covered " \
"(need %f%%)" % (grid_covered_ratio * 100,
grid_coverage * 100)
raise RuntimeError(msg)
LOG.debug("EWA resampling found %f%% of the grid covered" %
(grid_covered_ratio * 100))
return res
def compute(self, data, cache_id=None, fill_value=0, weight_count=10000,
weight_min=0.01, weight_distance_max=1.0,
weight_delta_max=1.0, weight_sum_min=-1.0,
maximum_weight_mode=False, grid_coverage=0, **kwargs):
"""Resample the data according to the precomputed X/Y coordinates."""
rows = self.cache["rows"]
cols = self.cache["cols"]
# if the data is scan based then check its metadata or the passed
# kwargs otherwise assume the entire input swath is one large
# "scanline"
rows_per_scan = kwargs.get('rows_per_scan',
data.attrs.get("rows_per_scan",
data.shape[0]))
if data.ndim == 3 and 'bands' in data.dims:
data_in = tuple(data.sel(bands=band).data
for band in data['bands'])
elif data.ndim == 2:
data_in = data.data
else:
raise ValueError("Unsupported data shape for EWA resampling.")
res = dask.delayed(self._call_fornav)(
cols, rows, self.target_geo_def, data_in,
grid_coverage=grid_coverage,
rows_per_scan=rows_per_scan, weight_count=weight_count,
weight_min=weight_min, weight_distance_max=weight_distance_max,
weight_delta_max=weight_delta_max, weight_sum_min=weight_sum_min,
maximum_weight_mode=maximum_weight_mode)
if isinstance(data_in, tuple):
new_shape = (len(data_in),) + self.target_geo_def.shape
else:
new_shape = self.target_geo_def.shape
data_arr = da.from_delayed(res, new_shape, data.dtype)
# from delayed creates one large chunk, break it up a bit if we can
data_arr = data_arr.rechunk([CHUNK_SIZE] * data_arr.ndim)
if data.ndim == 3 and data.dims[0] == 'bands':
dims = ('bands', 'y', 'x')
elif data.ndim == 2:
dims = ('y', 'x')
else:
dims = data.dims
res = xr.DataArray(data_arr, dims=dims, attrs=data.attrs.copy())
return update_resampled_coords(data, res, self.target_geo_def)
class BilinearResampler(BaseResampler):
"""Resample using bilinear interpolation.
This resampler implements on-disk caching when the `cache_dir` argument
is provided to the `resample` method. This should provide significant
performance improvements on consecutive resampling of geostationary data.
Args:
cache_dir (str): Long term storage directory for intermediate
results.
radius_of_influence (float): Search radius cut off distance in meters
epsilon (float): Allowed uncertainty in meters. Increasing uncertainty
reduces execution time.
reduce_data (bool): Reduce the input data to (roughly) match the
target area.
"""
def __init__(self, source_geo_def, target_geo_def):
"""Init BilinearResampler."""
super(BilinearResampler, self).__init__(source_geo_def, target_geo_def)
self.resampler = None
def precompute(self, mask=None, radius_of_influence=50000, epsilon=0,
reduce_data=True, cache_dir=False, **kwargs):
"""Create bilinear coefficients and store them for later use."""
try:
from pyresample.bilinear import XArrayBilinearResampler
except ImportError:
from pyresample.bilinear import XArrayResamplerBilinear as XArrayBilinearResampler
del kwargs
del mask
if self.resampler is None:
kwargs = dict(source_geo_def=self.source_geo_def,
target_geo_def=self.target_geo_def,
radius_of_influence=radius_of_influence,
neighbours=32,
epsilon=epsilon)
self.resampler = XArrayBilinearResampler(**kwargs)
try:
self.load_bil_info(cache_dir, **kwargs)
LOG.debug("Loaded bilinear parameters")
except IOError:
LOG.debug("Computing bilinear parameters")
self.resampler.get_bil_info()
LOG.debug("Saving bilinear parameters.")
self.save_bil_info(cache_dir, **kwargs)
def load_bil_info(self, cache_dir, **kwargs):
"""Load bilinear resampling info from cache directory."""
if cache_dir:
filename = self._create_cache_filename(cache_dir,
prefix='bil_lut-',
**kwargs)
try:
self.resampler.load_resampling_info(filename)
except AttributeError:
warnings.warn(
"Bilinear resampler can't handle caching, "
"please upgrade Pyresample to 0.17.0 or newer.",
stacklevel=2
)
raise IOError
else:
raise IOError
def save_bil_info(self, cache_dir, **kwargs):
"""Save bilinear resampling info to cache directory."""
if cache_dir:
filename = self._create_cache_filename(cache_dir,
prefix='bil_lut-',
**kwargs)
# There are some old caches, move them out of the way
if os.path.exists(filename):
_move_existing_caches(cache_dir, filename)
LOG.info('Saving BIL neighbour info to %s', filename)
try:
self.resampler.save_resampling_info(filename)
except AttributeError:
warnings.warn(
"Bilinear resampler can't handle caching, "
"please upgrade Pyresample to 0.17.0 or newer.",
stacklevel=2
)
def compute(self, data, fill_value=None, **kwargs):
"""Resample the given data using bilinear interpolation."""
del kwargs
if fill_value is None:
fill_value = data.attrs.get('_FillValue')
target_shape = self.target_geo_def.shape
res = self.resampler.get_sample_from_bil_info(data,
fill_value=fill_value,
output_shape=target_shape)
return update_resampled_coords(data, res, self.target_geo_def)
def _move_existing_caches(cache_dir, filename):
"""Move existing cache files out of the way."""
import os
import shutil
old_cache_dir = os.path.join(cache_dir, 'moved_by_satpy')
try:
os.makedirs(old_cache_dir)
except FileExistsError:
pass
try:
shutil.move(filename, old_cache_dir)
except shutil.Error:
os.remove(os.path.join(old_cache_dir,
os.path.basename(filename)))
shutil.move(filename, old_cache_dir)
LOG.warning("Old cache file was moved to %s", old_cache_dir)
def _mean(data, y_size, x_size):
rows, cols = data.shape
new_shape = (int(rows / y_size), int(y_size),
int(cols / x_size), int(x_size))
data_mean = np.nanmean(data.reshape(new_shape), axis=(1, 3))
return data_mean
def _repeat_by_factor(data, block_info=None):
if block_info is None:
return data
out_shape = block_info[None]['chunk-shape']
out_data = data
for axis, axis_size in enumerate(out_shape):
in_size = data.shape[axis]
out_data = np.repeat(out_data, int(axis_size / in_size), axis=axis)
return out_data
class NativeResampler(BaseResampler):
"""Expand or reduce input datasets to be the same shape.
If data is higher resolution (more pixels) than the destination area
then data is averaged to match the destination resolution.
If data is lower resolution (less pixels) than the destination area
then data is repeated to match the destination resolution.
This resampler does not perform any caching or masking due to the
simplicity of the operations.
"""