-
Notifications
You must be signed in to change notification settings - Fork 66
/
nansat.py
executable file
·2294 lines (1966 loc) · 87.9 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
# Name: nansat.py
# Purpose: Container of Nansat class
# Authors: Asuka Yamakawa, Anton Korosov, Knut-Frode Dagestad,
# Morten W. Hansen, Alexander Myasoyedov,
# Dmitry Petrenko, Evgeny Morozov
# Created: 29.06.2011
# Copyright: (c) NERSC 2011 - 2013
# 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
import os
import glob
import sys
import tempfile
import datetime
import dateutil.parser
import pkgutil
import warnings
try:
from collections import OrderedDict
except ImportError:
from ordereddict import OrderedDict
import scipy
from scipy.io.netcdf import netcdf_file
import numpy as np
import matplotlib
from matplotlib import cm
import matplotlib.pyplot as plt
from nansat.nsr import NSR
from nansat.domain import Domain
from nansat.figure import Figure
from nansat.vrt import VRT
from nansat.nansatshape import Nansatshape
from nansat.tools import add_logger, gdal
from nansat.tools import OptionError, WrongMapperError, NansatReadError, GDALError
from nansat.node import Node
from nansat.pointbrowser import PointBrowser
# container for all mappers
nansatMappers = None
def test_openable(fname):
try:
f = open(fname,'r')
except IOError:
raise
f.close()
class Nansat(Domain):
'''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.
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)
'''
def __init__(self, fileName='', mapperName='', domain=None,
array=None, parameters=None, logLevel=30, **kwargs):
'''Create Nansat object
if <fileName> is given:
Open GDAL dataset,
Read metadata,
Generate GDAL VRT file with mapping of variables in memory
Create logger
Create Nansat object for perfroming high-level operations
if <domain> and <array> are given:
Create VRT object from data in <array>
Add geolocation from <domain>
Parameters
-----------
fileName : string
location of the file
mapperName : string, optional
name of the mapper from nansat/mappers dir. E.g.
'ASAR', 'hirlam', 'merisL1', 'merisL2', etc.
domain : Domain object
Geo-reference of a new raster
array : numpy array
Firts band of a new raster
parameters : dictionary
Metadata for the 1st band of a new raster,e.g. name, wkv, units,...
logLevel : int, optional, default: logging.DEBUG (30)
Level of logging. See: http://docs.python.org/howto/logging.html
kwargs : additional arguments for mappers
Creates
--------
self.mapper : str
name of the used mapper
self.fileName : file name
set file name given by the argument
self.vrt : VRT object
Wrapper around VRT file and GDAL dataset with satellite raster data
self.logger : logging.Logger
logger for output debugging info
self.name : string
name of object (for writing KML)
Examples
--------
n = Nansat(filename)
# opens file for reading. Opening is lazy - no data is read at this
# point, only metadata that describes the dataset and bands
n = Nansat(domain=d)
# create an empty Nansat object. <d> is the Domain object which
# describes the grid (projection, resolution and extent)
n = Nansat(domain=d, array=a, parameters=p)
# create a Nansat object in memory with one band from input array <a>.
# <p> is a dictionary with metadata for the band
a = n[1]
# fetch data from Nansat object from the first band
a = n['band_name']
# fetch data from the band which has name 'band_name'
'''
# check the arguments
if fileName == '' and domain is None:
raise OptionError('Either fileName or domain is required.')
# create logger
self.logger = add_logger('Nansat', logLevel)
# empty dict of VRTs with added bands
self.addedBands = {}
# 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)
# create self.vrt from a file using mapper or...
if fileName != '':
# Make original VRT object with mapping of variables
self.vrt = self._get_mapper(mapperName, **kwargs)
# ...create using array, domain, and parameters
else:
# Set current VRT object
self.vrt = VRT(gdalDataset=domain.vrt.dataset)
self.domain = domain
self.mapper = ''
if array is not None:
# add a band from array
self.add_band(array=array, parameters=parameters)
self.logger.debug('Object created from %s ' % self.fileName)
def __getitem__(self, bandID):
''' Returns the band as a NumPy array, by overloading []
Parameters
-----------
bandID : int or str
If int, array from band with number <bandID> is returned
If string, array from band with metadata 'name' equal to
<bandID> is returned
Returns
--------
self.get_GDALRasterBand(bandID).ReadAsArray() : NumPy array
'''
# get band
band = self.get_GDALRasterBand(bandID)
# get expression from metadata
expression = band.GetMetadata().get('expression', '')
# get data
bandData = band.ReadAsArray()
if bandData is None:
raise GDALError('Cannot read array from band %s' % str(bandID))
# execute expression if any
if expression != '':
bandData = eval(expression)
# Set invalid and missing data to np.nan
if '_FillValue' in band.GetMetadata():
fillValue = float(band.GetMetadata()['_FillValue'])
try:
bandData[bandData == fillValue] = np.nan
except ValueError:
self.logger.info('Cannot replace _FillValue values '
'with np.NAN in %s!' % bandID)
try:
bandData[np.isinf(bandData)] = np.nan
except ValueError:
self.logger.info('Cannot replace inf values with np.NAN!')
return bandData
def __repr__(self):
'''Creates string with basic info about the Nansat object'''
outString = '-' * 40 + '\n'
outString += self.fileName + '\n'
outString += '-' * 40 + '\n'
outString += 'Mapper: ' + self.mapper + '\n'
outString += '-' * 40 + '\n'
outString += self.list_bands(False)
outString += '-' * 40 + '\n'
outString += Domain.__repr__(self)
return outString
def add_band(self, array, parameters=None, nomem=False):
'''Add band from the array to self.vrt
Create VRT object which contains VRT and RAW binary file and append it
to self.vrt.bandVRTs
Parameters
-----------
array : ndarray
band data
parameters : dictionary
band metadata: wkv, name, etc. (or for several bands)
nomem : boolean, saves the vrt to a tempfile if nomem is True
Modifies
---------
Creates VRT object with VRT-file and RAW-file
Adds band to the self.vrt
Examples
--------
n.add_band(a, p)
# add new band from numpy array <a> with metadata <p> in memory
# Shape of a should be equal to the shape of <n>
n.add_band(a, p, nomem=True)
# add new band from an array <a> with metadata <p> but keep it
# temporarli on disk intead of memory
'''
self.add_bands([array], [parameters], nomem)
def add_bands(self, arrays, parameters=None, nomem=False):
'''Add band from the array to self.vrt
Create VRT object which contains VRT and RAW binary file and append it
to self.vrt.bandVRTs
Parameters
-----------
array : ndarray or list
band data (or data for several bands)
parameters : dictionary or list
band metadata: wkv, name, etc. (or for several bands)
nomem : boolean, saves the vrt to a tempfile if nomem is True
Modifies
---------
Creates VRT object with VRT-file and RAW-file
Adds band to the self.vrt
Examples
--------
n.add_bands([a1, a2], [p1, p2])
# add two new bands from numpy arrays <a1> and <a2> with metadata in
# <p1> and <p2>
'''
# create VRTs from arrays
bandVRTs = [VRT(array=array, nomem=nomem) for array in arrays]
self.vrt = self.vrt.get_super_vrt()
# add the array band into self.vrt and get bandName
for bi, bandVRT in enumerate(bandVRTs):
params = parameters[bi]
if params is None:
params = {}
bandName = self.vrt._create_band(
{'SourceFilename': bandVRT.fileName,
'SourceBand': 1},
params)
self.vrt.bandVRTs[bandName] = bandVRT
self.vrt.dataset.FlushCache() # required after adding bands
def bands(self):
''' Make a dictionary with all metadata from all bands
Returns
--------
b : dictionary
key = N, value = dict with all band metadata
'''
b = {}
for iBand in range(self.vrt.dataset.RasterCount):
b[iBand + 1] = self.get_metadata(bandID=iBand + 1)
return b
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
'''
bandExists = False
for b in self.bands():
if self.bands()[b]['name'] == band or \
( self.bands()[b].has_key('standard_name') and
self.bands()[b]['standard_name'] == band):
bandExists = True
return bandExists
def export(self, fileName, bands=None, rmMetadata=[], addGeolocArray=True,
addGCPs=True, driver='netCDF', bottomup=False, options=None):
'''Export Nansat object into netCDF or GTiff file
Parameters
-----------
fileName : str
output file name
bands: list (default=None)
Specify band numbers to export.
If None, all bands are exported.
rmMetadata : list
metadata names for removal before export.
e.g. ['name', 'colormap', 'source', 'sourceBands']
addGeolocArray : bool
add geolocation array datasets to exported file?
addGCPs : bool
add GCPs to exported file?
driver : str
Name of GDAL driver (format)
bottomup : bool
False: Write swath-projected data with rows and columns organized
as in the original product.
True: Use the default behaviour of GDAL, which is to flip the rows
options : str or list
GDAL export options in format of: 'OPT=VAL', or
['OPT1=VAL1', 'OP2='VAL2']
See also http://www.gdal.org/frmt_netcdf.html
Modifies
---------
Create a netCDF file
!! NB
------
If number of bands is more than one,
serial numbers are added at the end of each band name.
It is possible to fix it by changing
line.4605 in GDAL/frmts/netcdf/netcdfdataset.cpp :
'if( nBands > 1 ) sprintf(szBandName,"%s%d",tmpMetadata,iBand);'
--> 'if( nBands > 1 ) sprintf(szBandName,"%s",tmpMetadata);'
CreateCopy fails in case the band name has special characters,
e.g. the slash in 'HH/VV'.
Examples
--------
n.export(netcdfile)
# export all the bands into a netDCF 3 file
n.export(driver='GTiff')
# export all bands into a GeoTiff
'''
# temporary VRT for exporting
exportVRT = self.vrt.copy()
exportVRT.real = []
exportVRT.imag = []
# delete unnecessary bands
if bands is not None:
srcBands = np.arange(self.vrt.dataset.RasterCount) + 1
dstBands = np.array(bands)
mask = np.in1d(srcBands, dstBands)
rmBands = srcBands[mask==False]
exportVRT.delete_bands(rmBands.tolist())
# Find complex data band
complexBands = []
node0 = Node.create(exportVRT.read_xml())
for iBand in node0.nodeList('VRTRasterBand'):
dataType = iBand.getAttribute('dataType')
if dataType[0] == 'C':
complexBands.append(int(iBand.getAttribute('band')))
# if data includes complex data,
# create two bands from real and imaginary data arrays
if len(complexBands) != 0:
for i in complexBands:
bandMetadataR = self.get_metadata(bandID=i)
bandMetadataR.pop('dataType')
try:
bandMetadataR.pop('PixelFunctionType')
except:
pass
# Copy metadata and modify 'name' for real and imag bands
bandMetadataI = bandMetadataR.copy()
bandMetadataR['name'] = bandMetadataR.pop('name') + '_real'
bandMetadataI['name'] = bandMetadataI.pop('name') + '_imag'
# Create bands from the real and imaginary numbers
exportVRT.real.append(VRT(array=self[i].real))
exportVRT.imag.append(VRT(array=self[i].imag))
metaDict = [{'src': {
'SourceFilename': exportVRT.real[-1].fileName,
'SourceBand': 1},
'dst': bandMetadataR},
{'src': {
'SourceFilename': exportVRT.imag[-1].fileName,
'SourceBand': 1},
'dst': bandMetadataI}]
exportVRT._create_bands(metaDict)
# delete the complex bands
exportVRT.delete_bands(complexBands)
# add bands with geolocation arrays to the VRT
if addGeolocArray and len(exportVRT.geolocationArray.d) > 0:
exportVRT._create_band(
{'SourceFilename': self.vrt.geolocationArray.d['X_DATASET'],
'SourceBand': int(self.vrt.geolocationArray.d['X_BAND'])},
{'wkv': 'longitude',
'name': 'GEOLOCATION_X_DATASET'})
exportVRT._create_band(
{'SourceFilename': self.vrt.geolocationArray.d['Y_DATASET'],
'SourceBand': int(self.vrt.geolocationArray.d['Y_BAND'])},
{'wkv': 'latitude',
'name': 'GEOLOCATION_Y_DATASET'})
gcps = exportVRT.dataset.GetGCPs()
if addGCPs and len(gcps) > 0:
# add GCPs in VRT metadata and remove geotransform
exportVRT._add_gcp_metadata(bottomup)
exportVRT._remove_geotransform()
# add projection metadata
srs = self.vrt.dataset.GetProjection()
exportVRT.dataset.SetMetadataItem('NANSAT_Projection',
srs.replace(',',
'|').replace('"', '&'))
# add GeoTransform metadata
geoTransformStr = str(self.vrt.dataset.GetGeoTransform()).replace(',',
'|')
exportVRT.dataset.SetMetadataItem('NANSAT_GeoTransform',
geoTransformStr)
# manage metadata for each band
for iBand in range(exportVRT.dataset.RasterCount):
band = exportVRT.dataset.GetRasterBand(iBand + 1)
bandMetadata = band.GetMetadata()
# set NETCDF_VARNAME
try:
bandMetadata['NETCDF_VARNAME'] = bandMetadata['name']
except:
self.logger.warning('Unable to set NETCDF_VARNAME for band %d'
% (iBand + 1))
# remove unwanted metadata from bands
for rmMeta in rmMetadata:
try:
bandMetadata.pop(rmMeta)
except:
self.logger.info('Unable to remove metadata'
'%s from band %d' % (rmMeta, iBand + 1))
band.SetMetadata(bandMetadata)
# remove unwanted global metadata
globMetadata = exportVRT.dataset.GetMetadata()
for rmMeta in rmMetadata:
try:
globMetadata.pop(rmMeta)
except:
self.logger.info('Global metadata %s not found' % rmMeta)
exportVRT.dataset.SetMetadata(globMetadata)
# if output filename is same as input one...
if self.fileName == fileName:
numOfBands = self.vrt.dataset.RasterCount
# create VRT from each band and add it
for iBand in range(numOfBands):
vrt = VRT(array=self[iBand + 1])
self.add_band(vrt=vrt)
metadata = self.get_metadata(bandID=iBand + 1)
self.set_metadata(key=metadata,
bandID=numOfBands + iBand + 1)
# remove source bands
self.vrt.delete_bands(range(1, numOfBands))
# get CreateCopy() options
if options is None:
options = []
if type(options) == str:
options = [options]
# set bottomup option
if bottomup:
options += ['WRITE_BOTTOMUP=NO']
else:
options += ['WRITE_BOTTOMUP=YES']
# Create an output file using GDAL
self.logger.debug('Exporting to %s using %s and %s...' % (fileName,
driver,
options))
dataset = gdal.GetDriverByName(driver).CreateCopy(fileName,
exportVRT.dataset,
options=options)
self.logger.debug('Export - OK!')
def export2thredds(self, fileName, bands, metadata=None,
maskName=None, rmMetadata=[],
time=None, createdTime=None):
''' Export data into a netCDF formatted for THREDDS server
Parameters
-----------
fileName : str
output file name
bands : dict
{'band_name': {'type' : '>i1',
'scale' : 0.1,
'offset' : 1000,
'metaKey1' : 'meta value 1',
'metaKey2' : 'meta value 2'}}
dictionary sets parameters for band creation
'type' - string representation of data type in the output band
'scale' - sets scale_factor and applies scaling
'offset' - sets 'scale_offset and applies offsetting
other entries (e.g. 'units': 'K') set other metadata
metadata : dict
Glbal metadata to add
maskName: string;
if data include a mask band: give the mask name.
Non-masked value is 64.
if None: no mask is added
rmMetadata : list
unwanted metadata names which will be removed
time : list with datetime objects
aqcuisition time of original data. That value will be in time dim
createdTime : datetime
date of creation. Will be in metadata 'created'
!! NB
------
Nansat object (self) has to be projected (with valid GeoTransform and
valid Spatial reference information) but not wth GCPs
Examples
--------
# create THREDDS formatted netcdf file with all bands and time variable
n.export2thredds(filename)
# export only the first band and add global metadata
n.export2thredds(filename, ['L_469'], {'description': 'example'})
# export several bands and modify type, scale and offset
bands = {'L_645' : {'type': '>i2', 'scale': 0.1, 'offset': 0},
'L_555' : {'type': '>i2', 'scale': 0.1, 'offset': 0}}
n.export2thredds(filename, bands)
'''
# raise error if self is not projected (has GCPs)
if len(self.vrt.dataset.GetGCPs()) > 0:
raise OptionError('Cannot export dataset with GCPS for THREDDS!')
# replace bands as list with bands as dict
if type(bands) is list:
bands = dict.fromkeys(bands, {})
# Create temporary empty Nansat object with self domain
data = Nansat(domain=self)
# get mask (if exist)
if maskName is not None:
mask = self[maskName]
# add required bands to data
dstBands = {}
srcBands = [self.bands()[b]['name'] for b in self.bands()]
for iband in bands:
# skip non exiting bands
if iband not in srcBands:
self.logger.error('%s is not found' % str(iband))
continue
array = self[iband]
# catch None band error
if array is None:
raise GDALError('%s is None' % str(iband))
# set type, scale and offset from input data or by default
dstBands[iband] = {}
dstBands[iband]['type'] = bands[iband].get('type',
array.dtype.str.replace('u', 'i'))
dstBands[iband]['scale'] = float(bands[iband].get('scale', 1.0))
dstBands[iband]['offset'] = float(bands[iband].get('offset', 0.0))
if '_FillValue' in bands[iband]:
dstBands[iband]['_FillValue'] = float(
bands[iband]['_FillValue'])
# mask values with np.nan
if maskName is not None and iband != maskName:
array[mask != 64] = np.nan
# add array to a temporary Nansat object
bandMetadata = self.get_metadata(bandID=iband)
# remove unwanted metadata from bands
for rmMeta in rmMetadata:
if rmMeta in bandMetadata.keys():
bandMetadata.pop(rmMeta)
data.add_band(array=array, parameters=bandMetadata)
self.logger.debug('Bands for export: %s' % str(dstBands))
# get corners of reprojected data
lonCrn, latCrn = data.get_corners()
# common global attributes:
if createdTime is None:
createdTime = (datetime.datetime.utcnow().
strftime('%Y-%m-%d %H:%M:%S UTC'))
globMetadata = {'institution': 'NERSC',
'source': 'satellite remote sensing',
'creation_date': createdTime,
'northernmost_latitude': np.float(max(latCrn)),
'southernmost_latitude': np.float(min(latCrn)),
'westernmost_longitude': np.float(min(lonCrn)),
'easternmost_longitude': np.float(max(lonCrn)),
'history': ' '}
#join or replace default by custom global metadata
if metadata is not None:
for metaKey in metadata:
globMetadata[metaKey] = metadata[metaKey]
# remove unwanted metadata from global metadata
for rmMeta in rmMetadata:
if rmMeta in globMetadata.keys():
globMetadata.pop(rmMeta)
# export temporary Nansat object to a temporary netCDF
fid, tmpName = tempfile.mkstemp(suffix='.nc')
data.export(tmpName)
# open files for input and output
ncI = netcdf_file(tmpName, 'r')
ncO = netcdf_file(fileName, 'w')
# collect info on dimention names
dimNames = []
gridMappingName = None
for ncIVarName in ncI.variables:
ncIVar = ncI.variables[ncIVarName]
dimNames += list(ncIVar.dimensions)
# get grid_mapping_name
if hasattr(ncIVar, 'grid_mapping_name'):
gridMappingName = ncIVar.grid_mapping_name
gridMappingVarName = ncIVarName
dimNames = list(set(dimNames))
# collect info on dimention shapes
dimShapes = {}
for dimName in dimNames:
dimVar = ncI.variables[dimName]
dimShapes[dimName] = dimVar.shape[0]
# create dimensions
for dimName in dimNames:
ncO.createDimension(dimName, dimShapes[dimName])
# add time dimention
ncO.createDimension('time', 1)
ncOVar = ncO.createVariable('time', '>f8', ('time', ))
ncOVar.calendar = 'standard'
ncOVar.long_name = 'time'
ncOVar.standard_name = 'time'
ncOVar.units = 'days since 1900-1-1 0:0:0 +0'
ncOVar.axis = 'T'
if time is None:
time = filter(None, self.get_time())
# create value of time variable
if len(time) > 0:
td = time[0] - datetime.datetime(1900, 1, 1)
days = td.days + (float(td.seconds) / 60.0 / 60.0 / 24.0)
# add date
ncOVar[:] = days
# recreate file
for ncIVarName in ncI.variables:
ncIVar = ncI.variables[ncIVarName]
self.logger.debug('Creating variable: %s' % ncIVarName)
if ncIVarName in ['x', 'y', 'lon', 'lat']:
# create simple x/y variables
ncOVar = ncO.createVariable(ncIVarName, '>f4',
ncIVar.dimensions)
elif ncIVarName == gridMappingVarName:
# create projection var
ncOVar = ncO.createVariable(ncIVarName, ncIVar.typecode(),
ncIVar.dimensions)
elif 'name' in ncIVar._attributes and ncIVar.name in dstBands:
# dont add time-axis to lon/lat grids
if ncIVar.name in ['lon', 'lat']:
dimensions = ncIVar.dimensions
else:
dimensions = ('time', ) + ncIVar.dimensions
ncOVar = ncO.createVariable(ncIVar.name,
dstBands[ncIVar.name]['type'],
dimensions)
# copy array from input data
data = np.array(ncIVar.data)
# copy rounded data from x/y
if ncIVarName in ['x', 'y']:
ncOVar[:] = np.floor(data).astype('>f4')
#add axis=X or axis=Y
ncOVar.axis = {'x': 'X', 'y': 'Y'}[ncIVarName]
for attrib in ncIVar._attributes:
if len(ncIVar._attributes[attrib]) > 0:
ncOVar._attributes[attrib] = ncIVar._attributes[attrib]
# copy data from lon/lat
if ncIVarName in ['lon', 'lat']:
ncOVar[:] = data.astype('>f4')
ncOVar._attributes = ncIVar._attributes
# copy projection data (only all attributes)
if ncIVarName == gridMappingVarName:
ncOVar._attributes = ncIVar._attributes
# copy data from variables in the list
if (len(ncIVar.dimensions) > 0 and
'name' in ncIVar._attributes and ncIVar.name in dstBands):
# add offset and scale attributes
scale = dstBands[ncIVar.name]['scale']
offset = dstBands[ncIVar.name]['offset']
if not (offset == 0.0 and scale == 1.0):
ncOVar._attributes['add_offset'] = offset
ncOVar._attributes['scale_factor'] = scale
data = (data - offset) / scale
# replace non-value by '_FillValue'
if (ncIVar.name in dstBands):
if '_FillValue' in dstBands[ncIVar.name].keys():
data[np.isnan(data)] = bands[ncIVar.name]['_FillValue']
ncOVar[:] = data.astype(dstBands[ncIVar.name]['type'])
# copy (some) attributes
for inAttrName in ncIVar._attributes:
if inAttrName not in ['dataType', 'SourceFilename',
'SourceBand', '_Unsigned',
'FillValue', 'time']:
ncOVar._attributes[inAttrName] = (
ncIVar._attributes[inAttrName])
# add custom attributes
if ncIVar.name in bands:
for newAttr in bands[ncIVar.name]:
if newAttr not in ['type', 'scale', 'offset']:
ncOVar._attributes[newAttr] = (
bands[ncIVar.name][newAttr])
# add grid_mapping info
if gridMappingName is not None:
ncOVar._attributes['grid_mapping'] = gridMappingName
# copy (some) global attributes
for globAttr in ncI._attributes:
if not(globAttr.strip().startswith('GDAL')):
ncO._attributes[globAttr] = ncI._attributes[globAttr]
# add common and custom global attributes
for globMeta in globMetadata:
ncO._attributes[globMeta] = globMetadata[globMeta]
# write output file
ncO.close()
# close original files
ncI.close()
# Delete the temprary netCDF file
fid = None
os.remove(tmpName)
return 0
def resize(self, factor=1, width=None, height=None,
pixelsize=None, eResampleAlg=-1):
'''Proportional resize of the dataset.
The dataset is resized as (xSize*factor, ySize*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.
eResampleAlg : int (GDALResampleAlg), optional
-1 : Average (default),
0 : NearestNeighbour
1 : Bilinear,
2 : Cubic,
3 : CubicSpline,
4 : Lancoz
Modifies
---------
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.
'''
# get current shape
rasterYSize = float(self.shape()[0])
rasterXSize = float(self.shape()[1])
# estimate factor if pixelsize is given
if pixelsize is not None:
deltaX, deltaY = self.get_pixelsize_meters()
factorX = deltaX / float(pixelsize)
factorY = deltaY / float(pixelsize)
factor = (factorX + factorY)/2
# estimate factor if width or height is given
if width is not None:
factor = float(width) / rasterXSize
if height is not None:
factor = float(height) / rasterYSize
# calculate new size
newRasterYSize = int(rasterYSize * factor)
newRasterXSize = int(rasterXSize * factor)
self.logger.info('New size/factor: (%f, %f)/%f' %
(newRasterXSize, newRasterYSize, factor))
if eResampleAlg <= 0:
self.vrt = self.vrt.get_subsampled_vrt(newRasterXSize,
newRasterYSize,
factor,
eResampleAlg)
else:
# update size and GeoTranform in XML of the warped VRT object
self.vrt = self.vrt.get_resized_vrt(newRasterXSize,
newRasterYSize,
eResampleAlg=eResampleAlg)
# 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 = 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, bandID=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
-----------
bandID : serial number or string, optional (default is 1)
if number - a band number of the band to fetch
if string bandID = {'name': bandID}
Returns
--------
GDAL RasterBand
Example
-------
b = n.get_GDALRasterBand(1)
b = n.get_GDALRasterBand('sigma0')
'''
# get band number
bandNumber = self._get_band_number(bandID)
# the GDAL RasterBand of the corresponding band is returned
return self.vrt.dataset.GetRasterBand(bandNumber)
def list_bands(self, doPrint=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
-----------
doPrint : boolean, optional, default=True
do print, otherwise it is returned as string
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 doPrint:
# print to screeen
print outString
else:
return outString
def reproject(self, dstDomain=None, eResampleAlg=0, blockSize=None,
WorkingDataType=None, tps=None, **kwargs):
''' Change projection of the object based on the given Domain
Create superVRT from self.vrt with AutoCreateWarpedVRT() using
projection from the dstDomain.
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 dstDomain is west of 0,
the object is shifted by 180 westwards.