-
Notifications
You must be signed in to change notification settings - Fork 52
/
Copy pathdata_handling.py
2063 lines (1666 loc) · 74.8 KB
/
data_handling.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: Data handling
# Purpose: Transforming netcdf to numpy arrays, checking mask file
#
# Author: PB
#
# Created: 13/07/2016
# Copyright: (c) PB 2016
# -------------------------------------------------------------------------
import os, glob
import calendar
#import numpy as np
from . import globals
from cwatm.management_modules.checks import *
from cwatm.management_modules.timestep import *
from cwatm.management_modules.replace_pcr import *
from cwatm.management_modules.messages import *
import difflib # to check the closest word in settingsfile, if an error occurs
import math
from cwatm.management_modules.dynamicModel import *
from netCDF4 import Dataset,num2date,date2num,date2index
#from netcdftime import utime
from osgeo import gdal
from osgeo import osr
from osgeo import gdalconst
import warnings
#import xarray as xr
#from concurrent.futures import ThreadPoolExecutor
#import asyncio
# -------------------------------------
"""
# Asychron reading of meteomaps to speed up - with the help of Jens de Bruijn
class AsyncReader:
def __init__(self, filepath, variable_name):
#def __init__(self,no):
#filepath = netcdf_file[no]
self.filepath = filepath
#variable_name = varnameGl[no]
#self.nf1 = Dataset(filepath,"r")
#self.ds = self.nf1.variables[variable_name]
self.ds = xr.open_dataset(filepath, chunks={"time": 1})[variable_name]
#self.ds = xr.open_dataarray(filepath)
# Adjust chunk size based on your data
self.executor = ThreadPoolExecutor(max_workers=1)
self.preloaded_data_future = None
self.current_index = -1 # Initialize to -1 to indicate no data loaded yet
#self.time_index = self.ds.time.values
#self.time_index = nf1.variables["time"][:]
try:
self.loop = asyncio.get_running_loop()
except RuntimeError:
self.loop = asyncio.new_event_loop()
def load(self, index,loc):
#a = self.ds[index + 1, loc[0]:loc[1], loc[2]:loc[3]]
#a = self.ds.values[index,loc[0]:loc[1],loc[2]:loc[3]]
#a = self.ds.isel(time=index).values
a = self.ds.isel(time=index, lat=slice(loc[0], loc[1]), lon=slice(loc[2],loc[3])).values
return a
def preload_next(self, index,loc):
# Preload the next timestep asynchronously
#if index + 1 < self.ds.time.size:
if index+1 < 14444111:
future = self.loop.run_in_executor(
self.executor, lambda: self.load(index + 1,loc) )
return future
return None
async def read_timestep_async(self, index, loc):
#assert index < self.ds.time.size, "Index out of bounds."
#assert index >= 0, "Index out of bounds."
# Check if the requested data is already preloaded, if so, just return that data
if index == self.current_index:
return self.current_data
# Check if the data for the next timestep is preloaded, if so, await for it to complete
if self.preloaded_data_future is not None and self.current_index + 1 == index:
data = await self.preloaded_data_future
# Load the requested data if not preloaded
else:
data = await self.loop.run_in_executor(
self.executor, lambda: self.load(index,loc))
# Initiate preloading the next timestep, do not await here, this returns a future
self.preloaded_data_future = self.preload_next(index,loc)
self.current_index = index
self.current_data = data
return data
def read_timestep(self, index, loc):
#index = self.get_index(date)
return self.loop.run_until_complete(self.read_timestep_async(index,loc))
def read_timestepxxx(self, index,loc):
#def read_timestep_not_async(self, index):
#index = self.get_index(date)
return self.load(index,loc)
#return self.ds[index,loc[0]:loc[1],loc[2]:loc[3]]
def close(self):
#self.nf1.close()
self.ds.close()
self.executor.shutdown()
"""
# --------------------------------------
def valuecell( coordx, coordstr, returnmap = True):
"""
to put a value into a raster map -> invert of cellvalue, map is converted into a numpy array first
:param coordx: x,y or lon/lat coordinate
:param coordstr: String of coordinates
:return: 1D array with new value
"""
coord = []
col = []
row = []
for xy in coordx:
try:
coord.append(float(xy))
except:
msg = "Error 101: Gauges in settings file: " + xy + " in " + coordstr + " is not a coordinate"
raise CWATMError(msg)
null = np.zeros((maskmapAttr['row'], maskmapAttr['col']))
null[null == 0] = -9999
for i in range(int(len(coord) / 2)):
col.append(int((coord[i * 2] - maskmapAttr['x']) * maskmapAttr['invcell']))
row.append(int((maskmapAttr['y'] - coord[i * 2 + 1]) * maskmapAttr['invcell']))
if col[i] >= 0 and row[i] >= 0 and col[i] < maskmapAttr['col'] and row[i] < maskmapAttr['row']:
null[row[i], col[i]] = i + 1
else:
x1 = maskmapAttr['x']
x2 = x1 + maskmapAttr['cell']* maskmapAttr['col']
y1 = maskmapAttr['y']
y2 = y1 - maskmapAttr['cell']* maskmapAttr['row']
box = "%5s %5.1f\n" %("",y1)
box += "%5s ---------\n" % ""
box += "%5s | |\n" % ""
box += "%5.1f | |%5.1f <-- Box of mask map\n" %(x1,x2)
box += "%5s | |\n" % ""
box += "%5s ---------\n" % ""
box += "%5s %5.1f\n" % ("", y2)
#print box
print("%2s %-17s %10s %8s" % ("No", "Name", "time[s]", "%"))
msg = "Error 102: Coordinates: x = " + str(coord[i * 2]) + ' y = ' + str(
coord[i * 2 + 1]) + " of gauge is outside mask map\n\n"
msg += box
msg +="\nPlease have a look at \"MaskMap\" or \"Gauges\""
raise CWATMError(msg)
if returnmap:
mapnp = compressArray(null).astype(np.int64)
return mapnp
else:
return col, row
def setmaskmapAttr(x,y,col,row,cell):
"""
Definition of cell size, coordinates of the meteo maps and maskmap
:param x: upper left corner x
:param y: upper left corner y
:param col: number of cols
:param row: number of rows
:param cell: cell size
:return: -
"""
invcell = round(1/cell,0)
# getgeotransform only delivers single precision!
if invcell == 0: invcell = 1/cell
cell = 1 / invcell
if (x-int(x)) != 0.:
if abs(x - int(x)) > 1e9:
x = 1/round(1/(x-int(x)),4) + int(x)
else: x = round(x,6)
if (y - int(y)) != 0.:
if abs(y - int(y)) > 1e9:
y = 1 / round(1 / (y - int(y)), 4) + int(y)
else: y = round(y,6)
# This is still not ok! Some rounding issues still appear sometimes
maskmapAttr['x'] = x
maskmapAttr['y'] = y
maskmapAttr['col'] = col
maskmapAttr['row'] = row
maskmapAttr['cell'] = cell
maskmapAttr['invcell'] = invcell
def loadsetclone(self,name):
"""
load the maskmap and set as clone
:param name: name of mask map, can be a file or - row col cellsize xupleft yupleft -
:return: new mask map
"""
filename = cbinding(name)
coord = filename.split()
if len(coord) == 2:
name = "Ldd"
if len(coord) == 5:
# changed order of x, y i- in setclone y is first in CWATM
# settings x is first
# setclone row col cellsize xupleft yupleft
# retancle: Number of Cols, Number of rows, cellsize, upper left corner X, upper left corner Y
mapnp = np.ones((int(coord[1]), int(coord[0])))
setmaskmapAttr(float(coord[3]),float(coord[4]), int(coord[0]),int(coord[1]),float(coord[2]))
#mapnp[mapnp == 0] = 1
#map = numpy2pcr(Boolean, mapnp, -9999)
elif len(coord) < 3:
filename = os.path.splitext(cbinding(name))[0] + '.nc'
try:
nf1 = Dataset(filename, 'r')
value = list(nf1.variables.items())[-1][0] # get the last variable name
x1 = list(nf1.variables.values())[0][0]
x2 = list(nf1.variables.values())[0][1]
xlast = list(nf1.variables.values())[0][-1]
#x1 = nf1.variables['lon'][0]
#x2 = nf1.variables['lon'][1]
#xlast = nf1.variables['lon'][-1]
#y1 = nf1.variables['lat'][0]
#ylast = nf1.variables['lat'][-1]
y1 = list(nf1.variables.values())[1][0]
ylast = list(nf1.variables.values())[1][-1]
# swap to make y1 the biggest number
if y1 < ylast: y1, ylast = ylast, y1
cellSize = np.abs(x2 - x1)
invcell = round(1/cellSize)
if invcell == 0: invcell = 1/cellSize
nrRows = int(0.5 + np.abs(ylast - y1) * invcell + 1)
nrCols = int(0.5 + np.abs(xlast - x1) * invcell + 1)
x = x1 - cellSize / 2
y = y1 + cellSize / 2
with warnings.catch_warnings():
warnings.simplefilter("ignore")
mapnp = np.array(nf1.variables[value][0:nrRows, 0:nrCols])
nf1.close()
setmaskmapAttr( x, y, nrCols, nrRows, cellSize)
flagmap = True
except:
# load geotiff
try:
filename = cbinding(name)
nf2 = gdal.Open(filename, gdalconst.GA_ReadOnly)
geotransform = nf2.GetGeoTransform()
geotrans.append(geotransform)
setmaskmapAttr( geotransform[0], geotransform[3], nf2.RasterXSize, nf2.RasterYSize, geotransform[1])
band = nf2.GetRasterBand(1)
#bandtype = gdal.GetDataTypeName(band.DataType)
mapnp = band.ReadAsArray(0, 0, nf2.RasterXSize, nf2.RasterYSize)
# 10 because that includes all valid LDD values [1-9]
mapnp[mapnp > 10] = 0
mapnp[mapnp < -10] = 0
flagmap = True
except:
raise CWATMFileError(filename,msg = "Error 201: File reading Error\n", sname=name)
if Flags['check']:
checkmap(name, filename, mapnp, flagmap, False,0)
else:
msg = "Error 103: Maskmap: " + filename + " is not a valid mask map nor valid coordinates nor valid point\n"
msg +="Or there is a whitespace or undefined character in Maskmap"
raise CWATMError(msg)
# put in the ldd map
# if there is no ldd at a cell, this cell should be excluded from modelling
maskldd = loadmap('Ldd', compress = False)
maskarea = np.bool_(mapnp)
mask = np.logical_not(np.logical_and(maskldd,maskarea))
# mask=np.isnan(mapnp)
# mask[mapnp==0] = True # all 0 become mask out
mapC = np.ma.compressed(np.ma.masked_array(mask,mask))
# Definition of compressed array and info how to blow it up again
maskinfo['mask']=mask
maskinfo['shape']=mask.shape
maskinfo['maskflat']=mask.ravel() # map to 1D not compresses
maskinfo['shapeflat']=maskinfo['maskflat'].shape #length of the 1D array
maskinfo['mapC']=mapC.shape # length of the compressed 1D array
maskinfo['maskall'] =np.ma.masked_all(maskinfo['shapeflat']) # empty map 1D but with mask
maskinfo['maskall'].mask = maskinfo['maskflat']
globals.inZero=np.zeros(maskinfo['mapC'])
if Flags['check']:
checkmap("Mask+Ldd", "", np.ma.masked_array(mask,mask), flagmap, True, mapC)
outpoints = 0
if len(coord) == 2:
outpoints = valuecell(coord, filename)
outpoints[outpoints < 0] = 0
print("Create catchment from point and river network")
mask2D, xleft, yup = self.routing_kinematic_module.catchment(outpoints)
mapC = maskfrompoint(mask2D, xleft, yup) + 1
area = np.sum(loadmap('CellArea')) * 1e-6
print("Number of cells in catchment: %6i = %7.0f km2" % (np.sum(mask2D), area))
# if the final results map should be cover up with some mask:
if "coverresult" in binding:
coverresult[0] = returnBool('coverresult')
if coverresult[0]:
cover = loadmap('covermap', compress=False, cut = False)
cover[cover > 1] = False
cover[cover == 1] = True
coverresult[1] = cover
return mapC
def maskfrompoint(mask2D, xleft, yup):
"""
load a static map either value or pc raster map or netcdf
:param mask2D: 2D array of new mask
:param xleft: left lon coordinate
:param yup: upper lat coordinate
:return: new mask map
"""
if xleft == -1:
msg = "Error 104: MaskMap point does not have a valid value in the river network (LDD)"
raise CWATMError(msg)
x = xleft * maskmapAttr['cell'] + maskmapAttr['x']
y = maskmapAttr['y'] - yup * maskmapAttr['cell']
maskmapAttr['x'] = x
maskmapAttr['y'] = y
maskmapAttr['col'] = mask2D.shape[1]
maskmapAttr['row'] = mask2D.shape[0]
mask = np.invert(np.bool_(mask2D))
mapC = np.ma.compressed(np.ma.masked_array(mask, mask))
# Definition of compressed array and info how to blow it up again
maskinfo['mask'] = mask
maskinfo['shape'] = mask.shape
maskinfo['maskflat'] = mask.ravel() # map to 1D not compresses
maskinfo['shapeflat'] = maskinfo['maskflat'].shape # length of the 1D array
maskinfo['mapC'] = mapC.shape # length of the compressed 1D array
maskinfo['maskall'] = np.ma.masked_all(maskinfo['shapeflat']) # empty map 1D but with mask
maskinfo['maskall'].mask = maskinfo['maskflat']
globals.inZero = np.zeros(maskinfo['mapC'])
return mapC
def loadmap(name, lddflag=False,compress = True, local = False, cut = True):
"""
load a static map either value or pc raster map or netcdf
:param name: name of map
:param lddflag: if True the map is used as a ldd map
:param compress: if True the return map will be compressed
:param local: if True the map is local and will be not cut
:param cut: if True the map will be not cut
:return: 1D numpy array of map
"""
value = cbinding(name)
filename = value
mapC = 0 # initializing to prevent warning in code inspection
try: # loading an integer or float but not a map
mapC = float(value)
flagmap = False
load = True
if Flags['check']:
checkmap(name, filename, mapC, False, False, 0)
except ValueError:
load = False
if not load: # read a netcdf (single one not a stack)
filename = os.path.splitext(value)[0] + '.nc'
# get mapextend of netcdf map and calculate the cutting
#cut0, cut1, cut2, cut3 = mapattrNetCDF(filename)
try:
nf1 = Dataset(filename, 'r')
cut0, cut1, cut2, cut3 = mapattrNetCDF(filename, check = False)
# load netcdf map but only the rectangle needed
#nf1 = Dataset(filename, 'r')
value = list(nf1.variables.items())[-1][0] # get the last variable name
if (nf1.variables[maskmapAttr['coordy']][0] - nf1.variables[maskmapAttr['coordy']][-1]) < 0:
msg = "Error 202: Latitude is in wrong order\n"
raise CWATMFileError(filename, msg)
if not timestepInit:
#with np.errstate(invalid='ignore'):
with warnings.catch_warnings():
warnings.simplefilter("ignore")
# in order to ignore some invalid value comments
if cut:
mapnp = nf1.variables[value][cut2:cut3, cut0:cut1].astype(np.float64)
else:
mapnp = nf1.variables[value][:]
else:
if 'time' in nf1.variables:
timestepI = Calendar(timestepInit[0])
if type(timestepI) is datetime.datetime:
timestepI = date2num(timestepI,nf1.variables['time'].units)
else: timestepI = int(timestepI) -1
if not(timestepI in nf1.variables['time'][:]):
msg = "Error 105 time step " + str(int(timestepI)+1)+" not stored in "+ filename
raise CWATMError(msg)
itime = np.where(nf1.variables['time'][:] == timestepI)[0][0]
if cut:
mapnp = nf1.variables[value][itime,cut2:cut3, cut0:cut1]
else:
mapnp = nf1.variables[value][itime][:]
else:
if cut:
mapnp = nf1.variables[value][cut2:cut3, cut0:cut1]
else:
mapnp = nf1.variables[value][:]
nf1.close()
except:
filename = cbinding(name)
try:
nf2 = gdal.Open(filename, gdalconst.GA_ReadOnly)
band = nf2.GetRasterBand(1)
mapnp = band.ReadAsArray(0, 0, nf2.RasterXSize, nf2.RasterYSize).astype(np.float64)
# if local no cut
if not local:
if cut:
cut0, cut1, cut2, cut3 = mapattrTiff(nf2)
mapnp = mapnp[cut2:cut3, cut0:cut1]
except:
msg = "Error 203: File does not exists"
raise CWATMFileError(filename,msg,sname=name)
try:
if any(maskinfo) and compress: mapnp.mask = maskinfo['mask']
except:
ii=0
if compress:
mapC = compressArray(mapnp,name=filename)
if Flags['check']:
checkmap(name, filename, mapnp, True, True, mapC)
else:
mapC = mapnp
if Flags['check']:
checkmap(name, filename, mapnp, True, False, 0)
return mapC
# -----------------------------------------------------------------------
# Compressing to 1-dimensional numpy array
# -----------------------------------------------------------------------
def compressArray(map, name="None", zeros = 0.):
"""
Compress 2D array with missing values to 1D array without missing values
:param map: in map
:param name: filename of the map
:param zeros: add zeros (default= 0) if values of map are to big or too small
:return: Compressed 1D array
"""
if map.shape != maskinfo['mask'].shape:
msg = "Error 105: " + name + " has a different shape than area or ldd \n"
raise CWATMError(msg)
mapnp1 = np.ma.masked_array(map, maskinfo['mask'])
mapC = np.ma.compressed(mapnp1)
# if fill: mapC[np.isnan(mapC)]=0
if name != "None":
if np.max(np.isnan(mapC)):
msg = "Error 106:" + name + " has less valid pixels than area or ldd \n"
raise CWATMError(msg)
# test if map has less valid pixel than area.map (or ldd)
# if a value is bigger or smaller than 1e20, -1e20 than the standard value is taken
mapC[mapC > 1.E20] = zeros
mapC[mapC < -1.E20] = zeros
return mapC
def decompress(map):
"""
Decompress 1D array without missing values to 2D array with missing values
:param map: numpy 1D array as input
:return: 2D array for displaying
"""
# dmap=np.ma.masked_all(maskinfo['shapeflat'], dtype=map.dtype)
dmap = maskinfo['maskall'].copy()
dmap[~maskinfo['maskflat']] = map[:]
dmap = dmap.reshape(maskinfo['shape'])
# check if integer map (like outlets, lakes etc
try:
checkint = str(map.dtype)
except:
checkint = "x"
if checkint == "int16" or checkint == "int32":
dmap[dmap.mask] = -9999
elif checkint == "int8":
dmap[dmap < 0] = 0
else:
dmap[dmap.mask] = -9999
return dmap
# -----------------------------------------------------------------------
# NETCDF
# -----------------------------------------------------------------------
def getmeta(key,varname,alternative):
"""
get the meta data information for the netcdf output from the global
variable metaNetcdfVar
:param key: key
:param varname: variable name e.g. self.var.Precipitation
:return: metadata information
"""
ret = alternative
if varname in metaNetcdfVar:
if key in metaNetcdfVar[varname]:
ret = metaNetcdfVar[varname][key]
return ret
def metaNetCDF():
"""
get the map metadata from precipitation netcdf maps
"""
try:
name = cbinding('PrecipitationMaps')
name1 = glob.glob(os.path.normpath(name))[0]
nf1 = Dataset(name1, 'r')
for var in nf1.variables:
metadataNCDF[var] = nf1.variables[var].__dict__
nf1.close()
except:
msg = "Error 204: Trying to get metadata from netcdf\n"
raise CWATMFileError(cbinding('PrecipitationMaps'),msg)
def readCoord(name):
"""
get the meta data information for the netcdf output from the global
variable metaNetcdfVar
:param name: name of the netcdf file
:return: latitude, longitude, cell size, inverse cell size
"""
namenc = os.path.splitext(name)[0] + '.nc'
try:
nf1 = Dataset(namenc, 'r')
nc = True
except:
nc = False
if nc:
lat, lon, cell, invcell, rows, cols = readCoordNetCDF(namenc)
else:
raster = gdal.Open(name)
rows = raster.RasterYSize
cols = raster.RasterXSize
gt = raster.GetGeoTransform()
cell = gt[1]
invcell = round(1.0 / cell, 0)
if invcell == 0: invcell = 1. / cell
# getgeotransform only delivers single precision!
cell = 1 / invcell
lon = gt[0]
lat = gt[3]
#lon = 1 / round(1 / (x1 - int(x1)), 4) + int(x1)
#lat = 1 / round(1 / (y1 - int(y1)), 4) + int(y1)
return lat, lon, cell, invcell, rows, cols
def readCoordNetCDF(name,check = True):
"""
reads the map attributes col, row etc from a netcdf map
:param name: name of the netcdf file
:param check: checking if netcdffile exists
:return: latitude, longitude, cell size, inverse cell size
:raises if no netcdf map can be found: :meth:`management_modules.messages.CWATMFileError`
"""
if check:
try:
nf1 = Dataset(name, 'r')
except:
msg = "Error 205: Checking netcdf map \n"
raise CWATMFileError(name,msg)
else:
# if subroutine is called already from inside a try command
nf1 = Dataset(name, 'r')
if not('coordx' in maskmapAttr.keys()):
if 'lon' in nf1.variables.keys():
maskmapAttr['coordx'] = 'lon'
maskmapAttr['coordy'] = 'lat'
else:
maskmapAttr['coordx'] = 'x'
maskmapAttr['coordy'] = 'y'
if 'X' in nf1.variables.keys():
maskmapAttr['coordx'] = 'X'
maskmapAttr['coordy'] = 'Y'
if 'x' in nf1.variables.keys():
maskmapAttr['coordx'] = 'x'
maskmapAttr['coordy'] = 'y'
rows = nf1.variables[maskmapAttr['coordy']].shape[0]
cols = nf1.variables[maskmapAttr['coordx']].shape[0]
lon0 = nf1.variables[maskmapAttr['coordx']][0]
lon1 = nf1.variables[maskmapAttr['coordx']][1]
lat0 = nf1.variables[maskmapAttr['coordy']][0]
latlast = nf1.variables[maskmapAttr['coordy']][-1]
nf1.close()
# swap to make lat0 the biggest number
if lat0 < latlast:
lat0, latlast = latlast, lat0
cell = round(np.abs(lon1 - lon0),8)
invcell = round(1.0 / cell, 0)
if invcell == 0: invcell = 1./cell
lon = round(lon0 - cell / 2,8)
lat = round(lat0 + cell / 2,8)
return lat,lon, cell,invcell,rows,cols
def readCalendar(name):
nf1 = Dataset(name, 'r')
dateVar['calendar'] = nf1.variables['time'].calendar
nf1.close()
def checkMeteo_Wordclim(meteodata, wordclimdata):
"""
reads the map attributes of meteo dataset and wordclima dataset
and compare if it has the same map extend
:param nmeteodata: name of the meteo netcdf file
:param wordlclimdata: cname of the wordlclim netcdf file
:return: True if meteo and wordclim has the same mapextend
:raises if map extend is different :meth:`management_modules.messages.CWATMFileError`
"""
try:
nf1 = Dataset(meteodata, 'r')
except:
msg = "Error 206: Checking netcdf map \n"
raise CWATMFileError(meteodata, msg)
if 'lon' in list(nf1.variables.keys()):
xy = ["lon", "lat"]
else:
xy = ["x", "y"]
lonM0 = nf1.variables[xy[0]][0]
lon1 = nf1.variables[xy[0]][1]
cellM = round(np.abs(lon1 - lonM0) / 2.,8)
lonM0 = round(lonM0 - cellM,8)
lonM1 = round(nf1.variables[xy[0]][-1] + cellM,8)
latM0 = nf1.variables[xy[1]][0]
latM1 = nf1.variables[xy[1]][-1]
nf1.close()
# swap to make lat0 the biggest number
if latM0 < latM1:
latM0, latM1 = latM1, latM0
latM0 = round(latM0 + cellM,8)
latM1 = round(latM1 - cellM,8)
# load Wordclima data
try:
nf1 = Dataset(wordclimdata, 'r')
except:
msg = "Error 207: Checking netcdf map \n"
raise CWATMFileError(wordclimdata, msg)
lonW0 = nf1.variables[xy[0]][0]
lon1 = nf1.variables[xy[0]][1]
cellW = round(np.abs(lon1 - lonW0) / 2.,8)
lonW0 = round(lonW0 - cellW,8)
lonW1 = round(nf1.variables[xy[0]][-1] + cellW,8)
latW0 = nf1.variables[xy[1]][0]
latW1 = nf1.variables[xy[1]][-1]
nf1.close()
# swap to make lat0 the biggest number
if latW0 < latW1:
latW0, latW1 = latW1, latW0
latW0 = round(latW0 + cellW,8)
latW1 = round(latW1 - cellW,8)
# calculate the controll variable
contr1 = (lonM0 + lonM1 + latM0 + latM1)
contr2 = (lonW0 + lonW1 + latW0 + latW1)
contr = abs(round(contr1 - contr2,5))
check = True
if contr > 0.00001:
#msg = "Data from meteo dataset and Wordclim dataset does not match"
#raise CWATMError(msg)
check = False
return check
def mapattrNetCDF(name, check=True):
"""
get the 4 corners of a netcdf map to cut the map
defines the rectangular of the mask map inside the netcdf map
calls function :meth:`management_modules.data_handling.readCoord`
:param name: name of the netcdf file
:param check: checking if netcdffile exists
:return: cut1,cut2,cut3,cut4
:raises if cell size is different: :meth:`management_modules.messages.CWATMError`
"""
lat, lon, cell, invcell, rows, cols = readCoord(name)
if maskmapAttr['invcell'] != invcell:
msg = "Error 107: Cell size different in maskmap: " + \
binding['MaskMap'] + " and: " + name
raise CWATMError(msg)
xx = maskmapAttr['x']
yy = maskmapAttr['y']
#cut0 = int(0.0001 + np.abs(xx - lon) * invcell) # argmin() ??
#cut2 = int(0.0001 + np.abs(yy - lat) * invcell) #
cut0 = int(np.abs(xx + maskmapAttr['cell']/2 - lon) * invcell) # argmin() ??
cut2 = int(np.abs(yy - maskmapAttr['cell']/2 - lat) * invcell) #
cut1 = cut0 + maskmapAttr['col']
cut3 = cut2 + maskmapAttr['row']
return cut0, cut1, cut2, cut3
def mapattrNetCDFMeteo(name, check = True):
"""
get the map attributes like col, row etc from a netcdf map
and define the rectangular of the mask map inside the netcdf map
calls function :meth:`management_modules.data_handling.readCoordNetCDF`
:param name: name of the netcdf file
:param check: checking if netcdffile exists
:return: cut0,cut1,cut2,cut3,cut4,cut5,cut6,cut7
"""
lat, lon, cell, invcell, rows, cols = readCoordNetCDF(name, check)
# x0,xend, y0,yend - borders of fine resolution map
lon0 = maskmapAttr['x']
lat0 = maskmapAttr['y']
lonend = lon0 + maskmapAttr['col'] / maskmapAttr['invcell']
latend = lat0 - maskmapAttr['row'] / maskmapAttr['invcell']
# cut for 0.5 deg map based on finer resolution
# lats = nc_simulated.variables['lat'][:]
# in_lat = discharge_location[1]
# lat_idx = geo_idx(in_lat, lats)
# geo_idx = (np.abs(dd_array - dd)).argmin()
# geo_idx(dd, dd_array):
cut0 = int(0.0001 + np.abs(lon0 - lon) * invcell)
cut2 = int(0.0001 + np.abs(lat0 - lat) * invcell)
# lon and lat of coarse meteo dataset
lonCoarse = (cut0 * cell) + lon
latCoarse = lat - (cut2 * cell)
cut4 = int(0.0001 + np.abs(lon0 - lonCoarse) * maskmapAttr['invcell'])
cut5 = cut4 + maskmapAttr['col']
cut6 = int(0.0001 + np.abs(lat0 - latCoarse) * maskmapAttr['invcell'])
cut7 = cut6 + maskmapAttr['row']
# now coarser cut of the coarse meteo dataset
cut1 = int(0.0001 + np.abs(lonend - lon) * invcell)
cut3 = int(0.0001 + np.abs(latend - lat) * invcell)
# test if fine cut is inside coarse cut
cellx = (cut1 - cut0) * maskmapAttr['reso_mask_meteo']
celly = (cut3 - cut2) * maskmapAttr['reso_mask_meteo']
if cellx < cut5:
cut1 += 1
if celly < cut7:
cut3 += 1
if maskmapAttr['coordy'] == 'lat':
if cut1 > (360 * invcell): cut1 = int(360 * invcell)
if cut3 > (180 * invcell): cut3 = int(180 * invcell)
return cut0, cut1, cut2, cut3, cut4, cut5, cut6, cut7
def mapattrTiff(nf2):
"""
map attributes of a geotiff file
:param nf2:
:return: cut0,cut1,cut2,cut3
"""
geotransform = nf2.GetGeoTransform()
x1 = geotransform[0]
y1 = geotransform[3]
#maskmapAttr['col'] = nf2.RasterXSize
#maskmapAttr['row'] = nf2.RasterYSize
cellSize = geotransform[1]
#invcell = round(1/cellSize,0)
if cellSize > 1:
invcell = 1 / cellSize
else:
invcell = round(1/cellSize,0)
# getgeotransform only delivers single precision!
cellSize = 1 / invcell
if (x1-int(x1)) != 0:
x1 = 1/round(1/(x1-int(x1)),4) + int(x1)
if (y1-int(y1)) != 0:
y1 = 1 / round(1 / (y1 - int(y1)), 4) + int(y1)
if maskmapAttr['invcell'] != invcell:
msg = "Error 108: Cell size different in maskmap: " + \
binding['MaskMap']
raise CWATMError(msg)
x = x1 - cellSize / 2
y = y1 + cellSize / 2
cut0 = int(0.01 + np.abs(maskmapAttr['x'] - x) * invcell)
cut2 = int(0.01 + np.abs(maskmapAttr['y'] - y) * invcell)
cut1 = cut0 + maskmapAttr['col']
cut3 = cut2 + maskmapAttr['row']
return cut0, cut1, cut2, cut3
def multinetdf(meteomaps, usebuffer,startcheck = 'dateBegin'):
"""
:param meteomaps: list of meteomaps to define start and end time
.param usebuffer: for certain map downscaling an additional cell around the maps is needed
:param startcheck: date of beginning simulation
:return
:raises if no map stack in meteo map folder: :meth:`management_modules.messages.CWATMFileError`
"""
end = dateVar['dateEnd']
no = 0
for maps in meteomaps:
name = cbinding(maps)
nameall = glob.glob(os.path.normpath(name))
if not nameall:
msg ="Error 208: File missing \n"
raise CWATMFileError(name,msg, sname=maps)
nameall.sort()
meteolist = {}
startfile = 0
for filename in nameall:
# --- Netcdf -------------
try:
nf1 = Dataset(filename, 'r')
except:
msg = "Error 209: Netcdf map stacks: " + filename +"\n"
raise CWATMFileError(filename, msg, sname=maps)
nctime = nf1.variables['time']
unitconv1 = ["DAYS", "HOUR", "MINU", "SECO"]
unitconv2 = [1, 24, 1440, 86400]
try:
unitconv3 = nctime.units[:4].upper()
datediv = unitconv2[unitconv1.index(unitconv3)]
except:
datediv = 1
datestart = num2date(int(round(nctime[:][0],0)), units=nctime.units,calendar=nctime.calendar)
# sometime daily records have a strange hour to start with -> it is changed to 0:00 to have the same record
datestart = datestart.replace(hour=0, minute=0)
datestartint = int(round(nctime[0].data.tolist(),0)) // datediv
datelen= len(nctime[:])
dateendint = datestartint + datelen - 1
dateend = num2date(dateendint, units=nctime.units, calendar=nctime.calendar)
#dateend = num2date(int(round(nctime[:][-1],0)), units=nctime.units, calendar=nctime.calendar)
#dateendint = int(round(nctime[:][-1].data.tolist(),0)) // datediv
dateend = dateend.replace(hour=0, minute=0)
#if dateVar['leapYear'] > 0:
startint = int(date2num(dateVar[startcheck],nctime.units,calendar=nctime.calendar))
start = num2date(startint, units=nctime.units, calendar=nctime.calendar)
startint = startint // datediv
endint = int(date2num(end, nctime.units, calendar=nctime.calendar))
endint = endint // datediv
#else:
# start = dateVar[startcheck]
value = list(nf1.variables.items())[-1][0] # get the last variable name
if value in ["X", "Y", "x", "y", "lon", "lat", "time"]:
for i in range(2, 5):
value = list(nf1.variables.items())[-i][0]
if not (value in ["X", "Y", "x", "y", "lon", "lat", "time"]): break
# check if mask = map size -> if yes do not cut the map
cutcheckmask = maskinfo['shape'][0] * maskinfo['shape'][1]
shapey = nf1.variables[value].shape[1]
shapex = nf1.variables[value].shape[2]
cutcheckmap = shapey * shapex
cutcheck = True
if cutcheckmask == cutcheckmap: cutcheck = False
# check if it is x or X
yy = maskmapAttr['coordy']
if yy == "y":
if "Y" in nf1.variables.keys():
yy = "Y"