-
Notifications
You must be signed in to change notification settings - Fork 66
/
nansat.py
1602 lines (1343 loc) · 58.6 KB
/
nansat.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
# Name: nansat.py
# Purpose: Container of Nansat class
# Authors: Anton Korosov, Knut-Frode Dagestad, Morten W. Hansen, Artem Moiseev,
# Asuka Yamakawa, Alexander Myasoyedov,
# Dmitry Petrenko, Evgeny Morozov
# Created: 29.06.2011
# Copyright: (c) NERSC 2011 - 2018
# Licence:
# This file is part of NANSAT.
# NANSAT 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, version 3 of the License.
# http://www.gnu.org/licenses/gpl-3.0.html
# 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.
from __future__ import absolute_import, print_function
import os
import glob
import sys
import tempfile
import datetime
import pkgutil
import warnings
from xml.sax import saxutils
import numpy as np
from numpy import nanmedian
from numpy.lib.recfunctions import append_fields
from netCDF4 import Dataset
from nansat.domain import Domain
from nansat.exporter import Exporter
from nansat.figure import Figure
from nansat.vrt import VRT
from nansat.tools import add_logger, gdal
from nansat.tools import parse_time, test_openable
from nansat.node import Node
from nansat.pointbrowser import PointBrowser
from nansat.exceptions import NansatGDALError, WrongMapperError, NansatReadError
import collections
if hasattr(collections, 'OrderedDict'):
from collections import OrderedDict
else:
from ordereddict import OrderedDict
# container for all mappers
nansatMappers = None
class Nansat(Domain, Exporter):
"""Container for geospatial data. Performs all high-level operations.
n = Nansat(filename) opens the file with satellite or model data for
reading, adds scientific metadata to bands, and prepares the data for
further handling.
Parameters
-----------
filename : str
uri of the input file or OpeNDAP datastream
mapper : str
name of the mapper from nansat/mappers dir. E.g.
'sentinel1_l1', 'asar', 'hirlam', 'meris_l1', 'meris_l2', etc.
log_level : int
Level of logging. See: http://docs.python.org/howto/logging.html
kwargs : additional arguments for mappers
Examples
--------
>>> n1 = Nansat(filename)
>>> n2 = Nansat(sentinel1_filename, mapper='sentinel1_l1')
>>> array1 = n1[1]
>>> array2 = n2['sigma0_HV']
Notes
-----
The instance of Nansat class (the object <n>) contains information
about geographical reference of the data (e.g raster size, pixel
resolution, type of projection, etc) and about bands with values of
geophysical variables (e.g. water leaving radiance, normalized radar
cross section, chlrophyll concentraion, etc). The object <n> has methods
for high-level operations with data. E.g.:
* reading data from file (Nansat.__getitem__);
* visualization (Nansat.write_figure);
* changing geographical reference (Nansat.reproject);
* exporting (Nansat.export)
* and much more...
Nansat inherits from Domain (container of geo-reference information)
Nansat uses instance of VRT (wraper around GDAL VRT-files)
Nansat uses instance of Figure (collection of methods for visualization)
"""
FILL_VALUE = 9.96921e+36
ALT_FILL_VALUE = -10000.
# instance attributes
logger = None
filename = None
name = None
path = None
vrt = None
mapper = None
@classmethod
def from_domain(cls, domain, array=None, parameters=None, log_level=30):
"""Create Nansat object from input Domain [and array with data]
Parameters
----------
domain : Domain
Defines spatial reference system and geographical extent.
array : numpy NDarray
Data for the first band. Shape must correspond to shape of <domain>
parameters : dict
Metadata for the first band. May contain 'name', 'wkv' and other keys.
log_level : int
Level of logging.
"""
n = cls.__new__(cls)
n._init_from_domain(domain, array, parameters, log_level)
return n
def __init__(self, filename='', mapper='', log_level=30, **kwargs):
"""Create Nansat object
Notes
-----
self.mapper : str
name of the used mapper
self.filename : file name
set file name given by the argument
self.vrt : VRT
Wrapper around VRT file and GDAL dataset with satellite raster data
self.logger : logging.Logger
logger for output debugging info
self.name : str
name of object (for writing KML)
self.path : str
path to input file
"""
if filename == '':
raise ValueError('Nansat is called without valid parameters! Use: Nansat(filename)')
self._init_empty(filename, log_level)
# Create VRT object with mapping of variables
self.vrt = self._get_mapper(mapper, **kwargs)
def __getitem__(self, band_id):
"""Returns the band as a NumPy array, by overloading []
Parameters
-----------
band_id : int or str
If int, array from band with number <band_id> is returned
If string, array from band with metadata 'name' equal to
<band_id> is returned
Returns
--------
a : NumPy array
"""
# get band
band = self.get_GDALRasterBand(band_id)
# get expression from metadata
expression = band.GetMetadata().get('expression', '')
# get data
band_data = band.ReadAsArray()
if band_data is None:
raise NansatGDALError('Cannot read array from band %s' % str(band_data))
# execute expression if any
if expression != '':
band_data = eval(expression)
all_float_flag = band_data.dtype.char in np.typecodes['AllFloat']
# Set invalid and missing data to np.nan (for floats only)
if '_FillValue' in band.GetMetadata() and all_float_flag:
band_data = self._fill_with_nan(band, band_data)
# replace infs with np.NAN
if np.size(np.where(np.isinf(band_data))) > 0:
band_data[np.isinf(band_data)] = np.nan
# erase out-of-swath pixels with np.Nan (if not integer)
if self.has_band('swathmask') and all_float_flag:
swathmask = self.get_GDALRasterBand('swathmask').ReadAsArray()
band_data[swathmask == 0] = np.nan
return band_data
def __repr__(self):
"""Creates string with basic info about the Nansat object"""
out_str = '{separator}{filename}{separator}Mapper: {mapper}{bands}{separator}{domain}'
return out_str.format(separator=self.OUTPUT_SEPARATOR, filename=self.filename,
bands=self.list_bands(False), mapper=self.mapper,
domain=Domain.__repr__(self))
def _init_empty(self, filename, log_level):
"""Init empty Nansat object
Parameters
----------
filename : str
Name of input file.
log_level : int
Level of logging verbosity.
Notes
--------
self.logger
adds logging.Logger
self.filename
adds full path to input file
self.name
adds name of file
self.path
adds path to file
"""
# create logger
self.logger = add_logger('Nansat', log_level)
# set input file name
self.filename = filename
# name, for compatibility with some Domain methods
self.name = os.path.basename(filename)
self.path = os.path.dirname(filename)
def _init_from_domain(self, domain, array=None, parameters=None, log_level=30):
"""Init Nansat object from input Domain and optionally array with band values
Parameters
----------
domain : Domain
Defines spatial reference system and geographical extent.
array : numpy NDarray
Data for the first band. Shape must correspond to shape of <domain>
parameters : dict
Metadata for the first band. May contain 'name', 'wkv' and other keys.
log_level : int
Level of logging.
"""
self._init_empty('', log_level)
self.vrt = VRT.from_gdal_dataset(domain.vrt.dataset, geolocation=domain.vrt.geolocation)
self.mapper = ''
if array is not None:
self.add_band(array=array, parameters=parameters)
def _fill_with_nan(self, band, band_data):
"""Fill input array with fill value taen from input band metadata"""
fill_value = float(band.GetMetadata()['_FillValue'])
band_data[band_data == fill_value] = np.nan
# quick hack to avoid problem with wrong _FillValue - see issue
# #123
if fill_value == self.FILL_VALUE:
band_data[band_data == self.ALT_FILL_VALUE] = np.nan
return band_data
def add_band(self, array, parameters=None, nomem=False):
"""Add band from numpy array with metadata.
Create VRT object which contains VRT and RAW binary file and append it
to self.vrt.band_vrts
Parameters
-----------
array : ndarray
new band data. Shape should be equal to shape
parameters : dict
band metadata: wkv, name, etc. (or for several bands)
nomem : bool
saves the vrt to a tempfile on disk?
Notes
-----
Creates VRT object with VRT-file and RAW-file. Adds band to the self.vrt.
Examples
--------
>>> n.add_band(array, {'name': 'new_data'}) # add new band and metadata, keep in memory
>>> n.add_band(array, nomem=True) # add new band, keep on disk
"""
self.add_bands([array], [parameters], nomem)
def add_bands(self, arrays, parameters=None, nomem=False):
"""Add bands from numpy arrays with metadata.
Create VRT object which contains VRT and RAW binary file and append it
to self.vrt.band_vrts
Parameters
-----------
arrays : list of ndarrays
new band data. Shape should be equal to shape
parameters : list of dict
band metadata: wkv, name, etc. (or for several bands)
nomem : bool
saves the vrt to a tempfile on disk?
Notes
-----
Creates VRT object with VRT-file and RAW-file. Adds band to the self.vrt.
Examples
--------
>>> n.add_bands([array1, array2]) # add new bands, keep in memory
"""
# replace empty parameters with list of empty dictionaries
if parameters is None:
parameters = [{}] * len(arrays)
self.vrt = self.vrt.get_super_vrt()
# create VRTs from arrays and generate band_metadata
band_metadata = []
for array, parameter in zip(arrays, parameters):
vrt = VRT.from_array(array, nomem=nomem)
band_metadata.append({
'src': {'SourceFilename': vrt.filename, 'SourceBand': 1},
'dst': parameter
})
self.vrt.band_vrts[vrt.filename] = vrt
band_name = self.vrt.create_bands(band_metadata)
def bands(self):
"""Make a dictionary with all metadata from all bands
Returns
--------
b : dictionary
key = N, value = dict with all band metadata
"""
band_metadata = {}
for band_num in range(self.vrt.dataset.RasterCount):
band_metadata[band_num + 1] = self.get_metadata(band_id=band_num + 1)
return band_metadata
def has_band(self, band):
"""Check if self has band with name <band>
Parameters
----------
band : str
name or standard_name of the band to check
Returns
-------
True/False if band exists or not
"""
for band_num in self.bands():
band_meta = self.bands()[band_num]
if band_meta['name'] == band:
return True
elif 'standard_name' in band_meta and band_meta['standard_name'] == band:
return True
def _get_resize_shape(self, factor, width, height, dst_pixel_size):
"""Estimate new shape either from factor or destination width/height or pixel size"""
src_shape = np.array(self.shape(), np.float)
# estimate factor if either width or height is given and factor is not given
if width is not None:
factor = width / src_shape[1]
if height is not None:
factor = height / src_shape[0]
# estimate factor if pixelsize is given
if dst_pixel_size is not None:
src_pixel_size = np.array(self.get_pixelsize_meters(), np.float)[::-1]
factor = (src_pixel_size / float(dst_pixel_size)).mean()
factor = float(factor)
dst_shape = np.floor(src_shape * factor)
self.logger.info('New shape: ({0}, {1})'.format(dst_shape[0], dst_shape[1]))
return factor, dst_shape
def resize(self, factor=None, width=None, height=None, pixelsize=None, resample_alg=-1):
"""Proportional resize of the dataset.
The dataset is resized as (x_size*factor, y_size*factor)
If desired width, height or pixelsize is specified,
the scaling factor is calculated accordingly.
If GCPs are given in a dataset, they are also rewritten.
Parameters
-----------
factor : float, optional, default=1
Scaling factor for width and height
> 1 means increasing domain size
< 1 means decreasing domain size
width : int, optional
Desired new width in pixels
height : int, optional
Desired new height in pixels
pixelsize : float, optional
Desired new pixelsize in meters (approximate).
A factor is calculated from ratio of the
current pixelsize to the desired pixelsize.
resample_alg : int (GDALResampleAlg), optional
-1 : Average (default),
0 : NearestNeighbour
1 : Bilinear,
2 : Cubic,
3 : CubicSpline,
4 : Lancoz
Notes
-----
self.vrt.dataset : VRT dataset of VRT object
raster size are modified to downscaled size.
If GCPs are given in the dataset, they are also overwritten.
"""
factor, dst_shape = self._get_resize_shape(factor, width, height, pixelsize)
if resample_alg <= 0:
self.vrt = self.vrt.get_subsampled_vrt(dst_shape[1], dst_shape[0], resample_alg)
else:
# update size and GeoTranform in XML of the warped VRT object
self.vrt = self.vrt.get_resized_vrt(dst_shape[1], dst_shape[0], resample_alg)
# resize gcps
gcps = self.vrt.vrt.dataset.GetGCPs()
if len(gcps) > 0:
gcpPro = self.vrt.vrt.dataset.GetGCPProjection()
for gcp in gcps:
gcp.GCPPixel *= factor
gcp.GCPLine *= factor
self.vrt.dataset.SetGCPs(gcps, gcpPro)
self.vrt._remove_geotransform()
else:
# change resultion in geotransform to keep spatial extent
geoTransform = list(self.vrt.vrt.dataset.GetGeoTransform())
geoTransform[1] = float(geoTransform[1])/factor
geoTransform[5] = float(geoTransform[5])/factor
geoTransform = list(map(float, geoTransform))
self.vrt.dataset.SetGeoTransform(geoTransform)
# set global metadata
subMetaData = self.vrt.vrt.dataset.GetMetadata()
subMetaData.pop('filename')
self.set_metadata(subMetaData)
return factor
def get_GDALRasterBand(self, band_id=1):
"""Get a GDALRasterBand of a given Nansat object
If str is given find corresponding band number
If int is given check if band with this number exists.
Get a GDALRasterBand from vrt.
Parameters
-----------
band_id : int or str
if int - a band number of the band to fetch
if str band_id = {'name': band_id}
Returns
--------
GDAL RasterBand
Example
-------
>>> b = n.get_GDALRasterBand(1)
>>> b = n.get_GDALRasterBand('sigma0')
"""
# get band number
bandNumber = self.get_band_number(band_id)
# the GDAL RasterBand of the corresponding band is returned
return self.vrt.dataset.GetRasterBand(bandNumber)
def list_bands(self, do_print=True):
"""Show band information of the given Nansat object
Show serial number, longName, name and all parameters
for each band in the metadata of the given Nansat object.
Parameters
-----------
do_print : boolean
print on screen?
Returns
--------
outString : String
formatted string with bands info
"""
# get dictionary of bands metadata
bands = self.bands()
outString = ''
for b in bands:
# print band number, name
outString += 'Band : %d %s\n' % (b, bands[b].get('name', ''))
# print band metadata
for i in bands[b]:
outString += ' %s: %s\n' % (i, bands[b][i])
if do_print:
# print to screeen
print(outString)
else:
return outString
def reproject(self, dst_domain=None, resample_alg=0,
block_size=None, tps=None, skip_gcps=1, addmask=True,
**kwargs):
"""
Change projection of the object based on the given Domain
Create superVRT from self.vrt with AutoCreateWarpedVRT() using
projection from the dst_domain.
Modify XML content of the warped vrt using the Domain parameters.
Generate warpedVRT and replace self.vrt with warpedVRT.
If current object spans from 0 to 360 and dst_domain is west of 0,
the object is shifted by 180 westwards.
Parameters
-----------
dst_domain : domain
destination Domain where projection and resolution are set
resample_alg : int (GDALResampleAlg)
0 : NearestNeighbour
1 : Bilinear
2 : Cubic,
3 : CubicSpline
4 : Lancoz
block_size : int
size of blocks for resampling. Large value decrease speed
but increase accuracy at the edge
tps : bool
Apply Thin Spline Transformation if source or destination has GCPs
Usage of TPS can also be triggered by setting self.vrt.tps=True
before calling to reproject.
This options has priority over self.vrt.tps
skip_gcps : int
Using TPS can be very slow if the number of GCPs are large.
If this parameter is given, only every [skip_gcp] GCP is used,
improving calculation time at the cost of accuracy.
If not given explicitly, 'skip_gcps' is fetched from the
metadata of self, or from dst_domain (as set by mapper or user).
[defaults to 1 if not specified, i.e. using all GCPs]
addmask : bool
If True, add band 'swathmask'. 1 - valid data, 0 no-data.
This band is used to replace no-data values with np.nan
Notes
-----
self.vrt : VRT object with dataset replaced to warpedVRT dataset
Integer data is returnd by integer. Round off to decimal place.
If you do not want to round off, convert the data types to
GDT_Float32, GDT_Float64, or GDT_CFloat32.
See Also
---------
http://www.gdal.org/gdalwarp.html
"""
# This is time consuming and therefore not done...:
#if not self.overlaps(dst_domain):
# raise ValueError('Source and destination domains do not overlap')
# if self spans from 0 to 360 and dst_domain is west of 0:
# shift self westwards by 180 degrees
# check span
srcCorners = self.get_corners()
if round(min(srcCorners[0])) == 0 and round(max(srcCorners[0])) == 360:
# check intersection of src and dst
dstCorners = dst_domain.get_corners()
if min(dstCorners[0]) < 0:
# shift
self.vrt = self.vrt.get_shifted_vrt(-180)
# get projection of destination dataset
dstSRS = dst_domain.vrt.dataset.GetProjection()
# get destination GCPs
dstGCPs = dst_domain.vrt.dataset.GetGCPs()
if len(dstGCPs) > 0:
# get projection of destination GCPs
dstSRS = dst_domain.vrt.dataset.GetGCPProjection()
x_size = dst_domain.vrt.dataset.RasterXSize
y_size = dst_domain.vrt.dataset.RasterYSize
geoTransform = dst_domain.vrt.dataset.GetGeoTransform()
# set trigger for using TPS
if tps is True:
self.vrt.tps = True
elif tps is False:
self.vrt.tps = False
# Reduce number of GCPs for faster reprojection
# when using TPS (if requested)
src_skip_gcps = self.vrt.dataset.GetMetadataItem('skip_gcps')
dst_skip_gcps = dst_domain.vrt.dataset.GetMetadataItem('skip_gcps')
kwargs['skip_gcps'] = skip_gcps # default (use all GCPs)
if dst_skip_gcps is not None: # ...or use setting from dst
kwargs['skip_gcps'] = int(dst_skip_gcps)
if src_skip_gcps is not None: # ...or use setting from src
kwargs['skip_gcps'] = int(src_skip_gcps)
# add band that masks valid values with 1 and nodata with 0
# after reproject
# TODD: REFACTOR: replace with VRT._add_swath_mask_band
if addmask:
self.vrt = self.vrt.get_super_vrt()
src = [{
'SourceFilename': self.vrt.vrt.filename,
'SourceBand': 1,
'DataType': gdal.GDT_Byte
}]
dst = {
'dataType': gdal.GDT_Byte,
'wkv': 'swath_binary_mask',
'PixelFunctionType': 'OnesPixelFunc',
}
self.vrt.create_band(src=src, dst=dst)
self.vrt.dataset.FlushCache()
# create Warped VRT
self.vrt = self.vrt.get_warped_vrt(dstSRS, x_size, y_size, geoTransform,
resample_alg=resample_alg,
dst_gcps=dstGCPs,
block_size=block_size, **kwargs)
# set global metadata from subVRT
subMetaData = self.vrt.vrt.dataset.GetMetadata()
subMetaData.pop('filename')
self.set_metadata(subMetaData)
def undo(self, steps=1):
"""Undo reproject, resize, add_band or crop of Nansat object
Restore the self.vrt from self.vrt.vrt
Parameters
-----------
steps : int
How many steps back to undo
Notes
------
Modifies self.vrt
"""
self.vrt = self.vrt.get_sub_vrt(steps)
def watermask(self, mod44path=None, dst_domain=None, **kwargs):
"""
Create numpy array with watermask (water=1, land=0)
250 meters resolution watermask from MODIS 44W Product:
http://www.glcf.umd.edu/data/watermask/
Watermask is stored as tiles in TIF(LZW) format and a VRT file
All files are stored in one directory.
A tarball with compressed TIF and VRT files should be additionally
downloaded from the Nansat documentation page:
http://nansat.readthedocs.io/en/latest/source/features.html#differentiating-between-land-and-water
The method :
Gets the directory either from input parameter or from environment
variable MOD44WPATH
Open Nansat object from the VRT file
Reprojects the watermask onto the current object using reproject()
or reproject_on_jcps()
Returns the reprojected Nansat object
Parameters
-----------
mod44path : string
path with MOD44W Products and a VRT file
dst_domain : Domain
destination domain other than self
tps : Bool
Use Thin Spline Transformation in reprojection of watermask?
See also Nansat.reproject()
skip_gcps : int
Factor to reduce the number of GCPs by and increase speed
See also Nansat.reproject()
Returns
--------
watermask : Nansat object with water mask in current projection
See Also
---------
http://www.glcf.umd.edu/data/watermask/
http://nansat.readthedocs.io/en/latest/source/features.html#differentiating-between-land-and-water
"""
mod44DataExist = True
# check if path is given in input param or in environment
if mod44path is None:
mod44path = os.getenv('MOD44WPATH')
if mod44path is None:
mod44DataExist = False
# check if VRT file exist
elif not os.path.exists(mod44path + '/MOD44W.vrt'):
mod44DataExist = False
self.logger.debug('MODPATH: %s' % mod44path)
if not mod44DataExist:
raise IOError('250 meters resolution watermask from MODIS '
'44W Product does not exist - see Nansat '
'documentation to get it (the path is % s)' % mod44path)
# MOD44W data does exist: open the VRT file in Nansat
watermask = Nansat(mod44path + '/MOD44W.vrt', mapper='MOD44W',
log_level=self.logger.level)
# reproject on self or given Domain
if dst_domain is None:
dst_domain = self
lon, lat = dst_domain.get_border()
watermask.crop_lonlat([lon.min(), lon.max()], [lat.min(), lat.max()])
watermask.reproject(dst_domain, addmask=False, **kwargs)
return watermask
def write_figure(self, filename='', bands=1, clim=None, addDate=False,
array_modfunc=None, **kwargs):
"""Save a raster band to a figure in graphical format.
Get numpy array from the band(s) and band information specified
either by given band number or band id.
-- If three bands are given, merge them and create PIL image.
-- If one band is given, create indexed image
Create Figure object and:
Adjust the array brightness and contrast using the given min/max or
histogram.
Apply logarithmic scaling of color tone.
Generate and append legend.
Save the PIL output image in PNG or any other graphical format.
If the filename extension is 'tif', the figure file is converted
to GeoTiff
Parameters
-----------
filename : str
Output file name. if one of extensions 'png', 'PNG', 'tif',
'TIF', 'bmp', 'BMP', 'jpg', 'JPG', 'jpeg', 'JPEG' is included,
specified file is created. otherwise, 'png' file is created.
bands : integer or string or list (elements are integer or string),
default = 1
the size of the list has to be 1 or 3.
if the size is 3, RGB image is created based on the three bands.
Then the first element is Red, the second is Green,
and the third is Blue.
clim : list with two elements or 'hist' to specify range of colormap
None (default) : min/max values are fetched from WKV,
fallback-'hist'
[min, max] : min and max are numbers, or
[[min, min, min], [max, max, max]]: three bands used
'hist' : a histogram is used to calculate min and max values
addDate : boolean
False (default) : no date will be aded to the caption
True : the first time of the object will be added to the caption
array_modfunc : None
None (default) : figure created using array in provided band
function : figure created using array modified by provided function
**kwargs : parameters for Figure().
Notes
---------
if filename is specified, creates image file
Returns
-------
Figure object
Example
--------
>>> n.write_figure('test.jpg') # write indexed image
>>> n.write_figure('test_rgb_hist.jpg', clim='hist', bands=[1, 2, 3]) # RGB image
>>> n.write_figure('r09_log3_leg.jpg', logarithm=True, legend=True,
gamma=3, titleString='Title', fontSize=30,
numOfTicks=15) # add legend
>>> n.write_figure(filename='transparent.png', bands=[3],
mask_array=wmArray,
mask_lut={0: [0,0,0]},
clim=[0,0.15], cmapName='gray',
transparency=[0,0,0]) # write transparent image
See also
--------
Figure()
http://www.scipy.org/Cookbook/Matplotlib/Show_colormaps
"""
# convert <bands> from integer, or string, or list of strings
# into list of integers
if isinstance(bands, list):
for i, band in enumerate(bands):
bands[i] = self.get_band_number(band)
else:
bands = [self.get_band_number(bands)]
# == create 3D ARRAY ==
array = None
for band in bands:
# get array from band and reshape to (1,height,width)
iArray = self[band]
if array_modfunc:
iArray = array_modfunc(iArray)
iArray = iArray.reshape(1, iArray.shape[0], iArray.shape[1])
# create new 3D array or append band
if array is None:
array = iArray
else:
array = np.append(array, iArray, axis=0)
# == CREATE FIGURE object and parse input parameters ==
fig = Figure(array, **kwargs)
array = None
# == PREPARE cmin/cmax ==
# check if cmin and cmax are given as the arguments
if 'cmin' in kwargs.keys() and 'cmax' in kwargs.keys():
clim = [kwargs['cmin'], kwargs['cmax']]
# try to get clim from WKV if it is not given as the argument
# if failed clim will be evaluated from histogram
if clim is None:
clim = [[], []]
for i, iBand in enumerate(bands):
try:
defValue = (self.vrt.dataset.GetRasterBand(iBand).
GetMetadataItem('minmax').split(' '))
except:
clim = 'hist'
break
clim[0].append(float(defValue[0]))
clim[1].append(float(defValue[1]))
# Estimate color min/max from histogram
if clim == 'hist':
clim = fig.clim_from_histogram(**kwargs)
# modify clim to the proper shape [[min], [max]]
# or [[min, min, min], [max, max, max]]
if (len(clim) == 2 and
((isinstance(clim[0], float)) or (isinstance(clim[0], int))) and
((isinstance(clim[1], float)) or (isinstance(clim[1], int)))):
clim = [[clim[0]], [clim[1]]]
# if the len(clim) is not same as len(bands), the 1st element is used.
for i in range(2):
if len(clim[i]) != len(bands):
clim[i] = [clim[i][0]] * len(bands)
self.logger.info('clim: %s ' % clim)
# == PREPARE caption ==
if 'caption' in kwargs:
caption = kwargs['caption']
else:
# get longName and units from vrt
band = self.get_GDALRasterBand(bands[0])
longName = band.GetMetadata().get('long_name', '')
units = band.GetMetadata().get('units', '')
# make caption from longname, units
caption = longName + ' [' + units + ']'
# add DATE to caption
if addDate:
caption += self.time_coverage_start.strftime(' %Y-%m-%d')
self.logger.info('caption: %s ' % caption)
# == PROCESS figure ==
fig.process(cmin=clim[0], cmax=clim[1], caption=caption)
# == finally SAVE to a image file ==
fig.save(filename, **kwargs)
# If tiff image, convert to GeoTiff
if filename[-3:] == 'tif':
self.vrt.copyproj(filename)
return fig
def write_geotiffimage(self, filename, band_id=1):
"""Writes an 8-bit GeoTiff image for a given band.
The colormap is fetched from the metadata item 'colormap'. Fallback colormap is 'gray'.
Color limits are fetched from the metadata item 'minmax'. If 'minmax' is not specified, min
and max of the raster data is used.
The method can be replaced by using nansat.write_figure(). However, write_figure uses PIL,
which does not allow Tiff compression. This gives much larger files.
Parameters
-----------
filename : str
band_id : int or str
"""
bandNo = self.get_band_number(band_id)
band = self.get_GDALRasterBand(band_id)
minmax = band.GetMetadataItem('minmax')
# Get min and max from band histogram if not given (from wkv)
if minmax is None:
(rmin, rmax) = band.ComputeRasterMinMax()
minmax = str(rmin) + ' ' + str(rmax)
bMin = float(minmax.split(' ')[0])
bMax = float(minmax.split(' ')[1])
# Make colormap from WKV information
cmap = np.vstack([np.arange(256.),
np.arange(256.),
np.arange(256.),
np.ones(256)*255]).T
colorTable = gdal.ColorTable()
for i in range(cmap.shape[0]):
colorEntry = (int(cmap[i, 0]), int(cmap[i, 1]),
int(cmap[i, 2]), int(cmap[i, 3]))
colorTable.SetColorEntry(i, colorEntry)
# Write Tiff image, with data scaled to values between 0 and 255
outDataset = gdal.GetDriverByName('Gtiff').Create(filename,
band.XSize,
band.YSize, 1,
gdal.GDT_Byte,
['COMPRESS=LZW'])
data = self.__getitem__(bandNo)
scaledData = ((data - bMin) / (bMax - bMin)) * 255
outDataset.GetRasterBand(1).WriteArray(scaledData)
outDataset.GetRasterBand(1).SetMetadata(band.GetMetadata())
try:
outDataset.GetRasterBand(1).SetColorTable(colorTable)
except:
# Happens after reprojection, a possible bug?
print('Could not set color table')
print(colorTable)
outDataset = None
self.vrt.copyproj(filename)
@property
def time_coverage_start(self):
return parse_time(self.get_metadata('time_coverage_start'))
@property
def time_coverage_end(self):
return parse_time(self.get_metadata('time_coverage_end'))
def get_metadata(self, key=None, band_id=None, unescape=True):
"""Get metadata from self.vrt.dataset
Parameters
----------
key : str
name of the metadata key. If not givem all metadata is returned
band_id : int or str
number or name of band to get metadata from.
If not given, global metadata is returned
unescape : bool
Replace '"', '&', '<' and '>' with these symbols " & < > ?
Returns
--------
* string with metadata if key is given and found
* dictionary with all metadata if key is not given
Raises
------
ValueError, if key is not found
"""
# get all metadata from dataset or from band
if band_id is None:
metadata = self.vrt.dataset.GetMetadata()
else:
metadata = self.get_GDALRasterBand(band_id).GetMetadata()
# remove escapes of special characters
if unescape:
for i in metadata:
metadata[i] = saxutils.unescape(metadata[i], {'"': '"'})
# get all metadata or from a key
if key is not None:
try: