forked from nasa-jpl/autoRIFT
-
Notifications
You must be signed in to change notification settings - Fork 2
/
testautoRIFT.py
executable file
·864 lines (690 loc) · 37.1 KB
/
testautoRIFT.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
#!/usr/bin/env python3
#~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
# Copyright 2019 California Institute of Technology. ALL RIGHTS RESERVED.
#
# Licensed under the Apache License, Version 2.0 (the "License");
# you may not use this file except in compliance with the License.
# You may obtain a copy of the License at
#
# http://www.apache.org/licenses/LICENSE-2.0
#
# Unless required by applicable law or agreed to in writing, software
# distributed under the License is distributed on an "AS IS" BASIS,
# WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
# See the License for the specific language governing permissions and
# limitations under the License.
#
# United States Government Sponsorship acknowledged. This software is subject to
# U.S. export control laws and regulations and has been classified as 'EAR99 NLR'
# (No [Export] License Required except when exporting to an embargoed country,
# end user, or in support of a prohibited end use). By downloading this software,
# the user agrees to comply with all applicable U.S. export laws and regulations.
# The user has the responsibility to obtain export licenses, or other export
# authority as may be required before exporting this software to any 'EAR99'
# embargoed foreign country or citizen of those countries.
#
# Author: Yang Lei
#~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
import pdb
from osgeo import gdal, osr
def runCmd(cmd):
import subprocess
out = subprocess.getoutput(cmd)
return out
def cmdLineParse():
'''
Command line parser.
'''
import argparse
parser = argparse.ArgumentParser(description='Output geo grid')
parser.add_argument('-m', '--input_m', dest='indir_m', type=str, required=True,
help='Input master image file name (in ISCE format and radar coordinates) or Input master image file name (in GeoTIFF format and Cartesian coordinates)')
parser.add_argument('-s', '--input_s', dest='indir_s', type=str, required=True,
help='Input slave image file name (in ISCE format and radar coordinates) or Input slave image file name (in GeoTIFF format and Cartesian coordinates)')
parser.add_argument('-g', '--input_g', dest='grid_location', type=str, required=False,
help='Input pixel indices file name')
parser.add_argument('-o', '--input_o', dest='init_offset', type=str, required=False,
help='Input search center offsets ("downstream" reach location) file name')
parser.add_argument('-sr', '--input_sr', dest='search_range', type=str, required=False,
help='Input search range file name')
parser.add_argument('-csmin', '--input_csmin', dest='chip_size_min', type=str, required=False,
help='Input chip size min file name')
parser.add_argument('-csmax', '--input_csmax', dest='chip_size_max', type=str, required=False,
help='Input chip size max file name')
parser.add_argument('-vx', '--input_vx', dest='offset2vx', type=str, required=False,
help='Input pixel offsets to vx conversion coefficients file name')
parser.add_argument('-vy', '--input_vy', dest='offset2vy', type=str, required=False,
help='Input pixel offsets to vy conversion coefficients file name')
parser.add_argument('-ssm', '--input_ssm', dest='stable_surface_mask', type=str, required=False,
help='Input stable surface mask file name')
parser.add_argument('-fo', '--flag_optical', dest='optical_flag', type=bool, required=False, default=0,
help='flag for reading optical data (e.g. Landsat): use 1 for on and 0 (default) for off')
parser.add_argument('-nc', '--sensor_flag_netCDF', dest='nc_sensor', type=str, required=False, default=None,
help='flag for packaging output formatted for Sentinel ("S") and Landsat ("L") dataset; default is None')
parser.add_argument('-urlflag', '--urlflag', dest='urlflag', type=int, required=False, default=0,
help='flag for reading and coregistering optical data (GeoTIFF images, e.g. Landsat): use 1 for url read and 0 (default) for local machine read')
parser.add_argument('-mpflag', '--mpflag', dest='mpflag', type=int, required=False, default=0,
help='number of threads for multiple threading (default is specified by 0, which uses the original single-core version and surpasses the multithreading routine)')
return parser.parse_args()
class Dummy(object):
pass
def loadProduct(filename):
'''
Load the product using Product Manager.
'''
import isce
import logging
from imageMath import IML
IMG = IML.mmapFromISCE(filename, logging)
img = IMG.bands[0]
# pdb.set_trace()
return img
def loadProductOptical(filename):
import numpy as np
'''
Load the product using Product Manager.
'''
ds = gdal.Open(filename)
# pdb.set_trace()
band = ds.GetRasterBand(1)
img = band.ReadAsArray()
img = img.astype(np.float32)
band=None
ds=None
return img
def loadProductOpticalURL(file_m, file_s):
import numpy as np
'''
Load the product using Product Manager.
'''
from geogrid import GeogridOptical
# from components.contrib.geo_autoRIFT.geogrid import GeogridOptical
obj = GeogridOptical()
x1a, y1a, xsize1, ysize1, x2a, y2a, xsize2, ysize2, trans = obj.coregister(file_m, file_s, 1)
DS1 = gdal.Open('/vsicurl/%s' %(file_m))
DS2 = gdal.Open('/vsicurl/%s' %(file_s))
I1 = DS1.ReadAsArray(xoff=x1a, yoff=y1a, xsize=xsize1, ysize=ysize1)
I2 = DS2.ReadAsArray(xoff=x2a, yoff=y2a, xsize=xsize2, ysize=ysize2)
I1 = I1.astype(np.float32)
I2 = I2.astype(np.float32)
DS1=None
DS2=None
return I1, I2
def runAutorift(I1, I2, xGrid, yGrid, Dx0, Dy0, SRx0, SRy0, CSMINx0, CSMINy0, CSMAXx0, CSMAXy0, noDataMask, optflag, nodata, mpflag):
'''
Wire and run geogrid.
'''
# import isce
# from components.contrib.geo_autoRIFT.autoRIFT import autoRIFT_ISCE
from autoRIFT import autoRIFT
import numpy as np
# import isceobj
import time
import subprocess
obj = autoRIFT()
# obj.configure()
# ########## uncomment if starting from preprocessed images
# I1 = I1.astype(np.uint8)
# I2 = I2.astype(np.uint8)
obj.MultiThread = mpflag
# take the amplitude only for the radar images
if optflag == 0:
I1 = np.abs(I1)
I2 = np.abs(I2)
obj.I1 = I1
obj.I2 = I2
# test with lena image (533 X 533)
# obj.ChipSizeMinX=16
# obj.ChipSizeMaxX=32
# obj.ChipSize0X=16
# obj.SkipSampleX=16
# obj.SkipSampleY=16
# test with Venus image (407 X 407)
# obj.ChipSizeMinX=8
# obj.ChipSizeMaxX=16
# obj.ChipSize0X=8
# obj.SkipSampleX=8
# obj.SkipSampleY=8
# create the grid if it does not exist
if xGrid is None:
m,n = obj.I1.shape
xGrid = np.arange(obj.SkipSampleX+10,n-obj.SkipSampleX,obj.SkipSampleX)
yGrid = np.arange(obj.SkipSampleY+10,m-obj.SkipSampleY,obj.SkipSampleY)
nd = xGrid.__len__()
md = yGrid.__len__()
obj.xGrid = np.int32(np.dot(np.ones((md,1)),np.reshape(xGrid,(1,xGrid.__len__()))))
obj.yGrid = np.int32(np.dot(np.reshape(yGrid,(yGrid.__len__(),1)),np.ones((1,nd))))
noDataMask = np.logical_not(obj.xGrid)
else:
obj.xGrid = xGrid
obj.yGrid = yGrid
# generate the nodata mask where offset searching will be skipped based on 1) imported nodata mask and/or 2) zero values in the image
for ii in range(obj.xGrid.shape[0]):
for jj in range(obj.xGrid.shape[1]):
if (obj.yGrid[ii,jj] != nodata)&(obj.xGrid[ii,jj] != nodata):
if (I1[obj.yGrid[ii,jj]-1,obj.xGrid[ii,jj]-1]==0)|(I2[obj.yGrid[ii,jj]-1,obj.xGrid[ii,jj]-1]==0):
noDataMask[ii,jj] = True
######### mask out nodata to skip the offset searching using the nodata mask (by setting SearchLimit to be 0)
if SRx0 is None:
# ########### uncomment to customize SearchLimit based on velocity distribution (i.e. Dx0 must not be None)
# obj.SearchLimitX = np.int32(4+(25-4)/(np.max(np.abs(Dx0[np.logical_not(noDataMask)]))-np.min(np.abs(Dx0[np.logical_not(noDataMask)])))*(np.abs(Dx0)-np.min(np.abs(Dx0[np.logical_not(noDataMask)]))))
# obj.SearchLimitY = 5
# ###########
obj.SearchLimitX = obj.SearchLimitX * np.logical_not(noDataMask)
obj.SearchLimitY = obj.SearchLimitY * np.logical_not(noDataMask)
else:
obj.SearchLimitX = SRx0
obj.SearchLimitY = SRy0
# ############ add buffer to search range
# obj.SearchLimitX[obj.SearchLimitX!=0] = obj.SearchLimitX[obj.SearchLimitX!=0] + 2
# obj.SearchLimitY[obj.SearchLimitY!=0] = obj.SearchLimitY[obj.SearchLimitY!=0] + 2
if CSMINx0 is not None:
obj.ChipSizeMaxX = CSMAXx0
obj.ChipSizeMinX = CSMINx0
chipsizex0 = float(str.split(subprocess.getoutput('fgrep "Smallest Allowable Chip Size in m:" testGeogrid.txt'))[-1])
try:
pixsizex = float(str.split(subprocess.getoutput('fgrep "Ground range pixel size:" testGeogrid.txt'))[-1])
except:
pixsizex = float(str.split(subprocess.getoutput('fgrep "X-direction pixel size:" testGeogrid.txt'))[-1])
obj.ChipSize0X = int(np.ceil(chipsizex0/pixsizex/4)*4)
# obj.ChipSize0X = np.min(CSMINx0[CSMINx0!=nodata])
RATIO_Y2X = CSMINy0/CSMINx0
obj.ScaleChipSizeY = np.mean(RATIO_Y2X[(CSMINx0!=nodata)&(CSMINy0!=nodata)])
# obj.ChipSizeMaxX = obj.ChipSizeMaxX / obj.ChipSizeMaxX * 544
# obj.ChipSizeMinX = obj.ChipSizeMinX / obj.ChipSizeMinX * 68
else:
if ((optflag == 1)&(xGrid is not None)):
# print('test')
obj.ChipSizeMaxX = 32
obj.ChipSizeMinX = 16
obj.ChipSize0X = 16
# create the downstream search offset if not provided as input
if Dx0 is not None:
obj.Dx0 = Dx0
obj.Dy0 = Dy0
else:
obj.Dx0 = obj.Dx0 * np.logical_not(noDataMask)
obj.Dy0 = obj.Dy0 * np.logical_not(noDataMask)
# replace the nodata value with zero
obj.xGrid[noDataMask] = 0
obj.yGrid[noDataMask] = 0
obj.Dx0[noDataMask] = 0
obj.Dy0[noDataMask] = 0
if SRx0 is not None:
obj.SearchLimitX[noDataMask] = 0
obj.SearchLimitY[noDataMask] = 0
if CSMINx0 is not None:
obj.ChipSizeMaxX[noDataMask] = 0
obj.ChipSizeMinX[noDataMask] = 0
# convert azimuth offset to vertical offset as used in autoRIFT convention
if optflag == 0:
obj.Dy0 = -1 * obj.Dy0
######## preprocessing
t1 = time.time()
print("Pre-process Start!!!")
# obj.zeroMask = 1
obj.preprocess_filt_hps()
# obj.I1 = np.abs(I1)
# obj.I2 = np.abs(I2)
print("Pre-process Done!!!")
print(time.time()-t1)
t1 = time.time()
# obj.DataType = 0
obj.uniform_data_type()
print("Uniform Data Type Done!!!")
print(time.time()-t1)
# pdb.set_trace()
# obj.sparseSearchSampleRate = 16
obj.OverSampleRatio = 64
# OverSampleRatio can be assigned as a scalar (such as the above line) or as a Python dictionary below for intellgient use (ChipSize-dependent).
# Here, four chip sizes are used: ChipSize0X*[1,2,4,8] and four OverSampleRatio are considered [16,32,64,128]. The intelligent selection of OverSampleRatio (as a function of chip size) was determined by analyzing various combinations of (OverSampleRatio and chip size) and comparing the resulting image quality and statistics with the reference scenario (where the largest OverSampleRatio of 128 and chip size of ChipSize0X*8 are considered).
# The selection for the optical data flag is based on Landsat-8 data over an inland region (thus stable and not moving much) of Greenland, while that for the radar flag (optflag = 0) is based on Sentinel-1 data over the same region of Greenland.
if CSMINx0 is not None:
if (optflag == 1):
obj.OverSampleRatio = {obj.ChipSize0X:16,obj.ChipSize0X*2:32,obj.ChipSize0X*4:64,obj.ChipSize0X*8:64}
else:
obj.OverSampleRatio = {obj.ChipSize0X:32,obj.ChipSize0X*2:64,obj.ChipSize0X*4:128,obj.ChipSize0X*8:128}
# ########## export preprocessed images to files; can be commented out if not debugging
#
# t1 = time.time()
#
# I1 = obj.I1
# I2 = obj.I2
#
# length,width = I1.shape
#
# filename1 = 'I1_uint8_hpsnew.off'
#
# slcFid = open(filename1, 'wb')
#
# for yy in range(length):
# data = I1[yy,:]
# data.astype(np.float32).tofile(slcFid)
#
# slcFid.close()
#
# img = isceobj.createOffsetImage()
# img.setFilename(filename1)
# img.setBands(1)
# img.setWidth(width)
# img.setLength(length)
# img.setAccessMode('READ')
# img.renderHdr()
#
#
# filename2 = 'I2_uint8_hpsnew.off'
#
# slcFid = open(filename2, 'wb')
#
# for yy in range(length):
# data = I2[yy,:]
# data.astype(np.float32).tofile(slcFid)
#
# slcFid.close()
#
# img = isceobj.createOffsetImage()
# img.setFilename(filename2)
# img.setBands(1)
# img.setWidth(width)
# img.setLength(length)
# img.setAccessMode('READ')
# img.renderHdr()
#
# print("output Done!!!")
# print(time.time()-t1)
########## run Autorift
t1 = time.time()
print("AutoRIFT Start!!!")
obj.runAutorift()
print("AutoRIFT Done!!!")
print(time.time()-t1)
import cv2
kernel = np.ones((3,3),np.uint8)
noDataMask = cv2.dilate(noDataMask.astype(np.uint8),kernel,iterations = 1)
noDataMask = noDataMask.astype(np.bool)
return obj.Dx, obj.Dy, obj.InterpMask, obj.ChipSizeX, obj.ScaleChipSizeY, obj.SearchLimitX, obj.SearchLimitY, obj.origSize, noDataMask
if __name__ == '__main__':
'''
Main driver.
'''
import numpy as np
import time
inps = cmdLineParse()
mpflag = inps.mpflag
urlflag = inps.urlflag
if inps.optical_flag == 1:
if urlflag == 1:
data_m, data_s = loadProductOpticalURL(inps.indir_m, inps.indir_s)
else:
data_m = loadProductOptical(inps.indir_m)
data_s = loadProductOptical(inps.indir_s)
# test with lena/Venus image
# import scipy.io as sio
# conts = sio.loadmat(inps.indir_m)
# data_m = conts['I']
# data_s = conts['I1']
else:
data_m = loadProduct(inps.indir_m)
data_s = loadProduct(inps.indir_s)
xGrid = None
yGrid = None
Dx0 = None
Dy0 = None
SRx0 = None
SRy0 = None
CSMINx0 = None
CSMINy0 = None
CSMAXx0 = None
CSMAXy0 = None
noDataMask = None
nodata = None
if inps.grid_location is not None:
ds = gdal.Open(inps.grid_location)
tran = ds.GetGeoTransform()
proj = ds.GetProjection()
srs = ds.GetSpatialRef()
band = ds.GetRasterBand(1)
nodata = band.GetNoDataValue()
xGrid = band.ReadAsArray()
noDataMask = (xGrid == nodata)
band = ds.GetRasterBand(2)
yGrid = band.ReadAsArray()
band=None
ds=None
if inps.init_offset is not None:
ds = gdal.Open(inps.init_offset)
band = ds.GetRasterBand(1)
Dx0 = band.ReadAsArray()
band = ds.GetRasterBand(2)
Dy0 = band.ReadAsArray()
band=None
ds=None
if inps.search_range is not None:
ds = gdal.Open(inps.search_range)
band = ds.GetRasterBand(1)
SRx0 = band.ReadAsArray()
band = ds.GetRasterBand(2)
SRy0 = band.ReadAsArray()
band=None
ds=None
if inps.chip_size_min is not None:
ds = gdal.Open(inps.chip_size_min)
band = ds.GetRasterBand(1)
CSMINx0 = band.ReadAsArray()
band = ds.GetRasterBand(2)
CSMINy0 = band.ReadAsArray()
band=None
ds=None
if inps.chip_size_max is not None:
ds = gdal.Open(inps.chip_size_max)
band = ds.GetRasterBand(1)
CSMAXx0 = band.ReadAsArray()
band = ds.GetRasterBand(2)
CSMAXy0 = band.ReadAsArray()
band=None
ds=None
if inps.stable_surface_mask is not None:
ds = gdal.Open(inps.stable_surface_mask)
band = ds.GetRasterBand(1)
SSM = band.ReadAsArray()
SSM = SSM.astype('bool')
band=None
ds=None
Dx, Dy, InterpMask, ChipSizeX, ScaleChipSizeY, SearchLimitX, SearchLimitY, origSize, noDataMask = runAutorift(data_m, data_s, xGrid, yGrid, Dx0, Dy0, SRx0, SRy0, CSMINx0, CSMINy0, CSMAXx0, CSMAXy0, noDataMask, inps.optical_flag, nodata, mpflag)
if inps.optical_flag == 0:
Dy = -Dy
DX = np.zeros(origSize,dtype=np.float32) * np.nan
DY = np.zeros(origSize,dtype=np.float32) * np.nan
INTERPMASK = np.zeros(origSize,dtype=np.float32)
CHIPSIZEX = np.zeros(origSize,dtype=np.float32)
SEARCHLIMITX = np.zeros(origSize,dtype=np.float32)
SEARCHLIMITY = np.zeros(origSize,dtype=np.float32)
DX[0:Dx.shape[0],0:Dx.shape[1]] = Dx
DY[0:Dy.shape[0],0:Dy.shape[1]] = Dy
INTERPMASK[0:InterpMask.shape[0],0:InterpMask.shape[1]] = InterpMask
CHIPSIZEX[0:ChipSizeX.shape[0],0:ChipSizeX.shape[1]] = ChipSizeX
SEARCHLIMITX[0:SearchLimitX.shape[0],0:SearchLimitX.shape[1]] = SearchLimitX
SEARCHLIMITY[0:SearchLimitY.shape[0],0:SearchLimitY.shape[1]] = SearchLimitY
DX[noDataMask] = np.nan
DY[noDataMask] = np.nan
INTERPMASK[noDataMask] = 0
CHIPSIZEX[noDataMask] = 0
SEARCHLIMITX[noDataMask] = 0
SEARCHLIMITY[noDataMask] = 0
if SSM is not None:
SSM[noDataMask] = False
import scipy.io as sio
sio.savemat('offset.mat',{'Dx':DX,'Dy':DY,'InterpMask':INTERPMASK,'ChipSizeX':CHIPSIZEX})
# ##################### Uncomment for debug mode
# sio.savemat('debug.mat',{'Dx':DX,'Dy':DY,'InterpMask':INTERPMASK,'ChipSizeX':CHIPSIZEX,'ScaleChipSizeY':ScaleChipSizeY,'SearchLimitX':SEARCHLIMITX,'SearchLimitY':SEARCHLIMITY})
# conts = sio.loadmat('debug.mat')
# DX = conts['Dx']
# DY = conts['Dy']
# INTERPMASK = conts['InterpMask']
# CHIPSIZEX = conts['ChipSizeX']
# ScaleChipSizeY = conts['ScaleChipSizeY']
# SEARCHLIMITX = conts['SearchLimitX']
# SEARCHLIMITY = conts['SearchLimitY']
# #####################
if inps.grid_location is not None:
t1 = time.time()
print("Write Outputs Start!!!")
# Create the GeoTiff
driver = gdal.GetDriverByName('GTiff')
outRaster = driver.Create("offset.tif", int(xGrid.shape[1]), int(xGrid.shape[0]), 4, gdal.GDT_Float32)
outRaster.SetGeoTransform(tran)
outRaster.SetProjection(proj)
outband = outRaster.GetRasterBand(1)
outband.WriteArray(DX)
outband.FlushCache()
outband = outRaster.GetRasterBand(2)
outband.WriteArray(DY)
outband.FlushCache()
outband = outRaster.GetRasterBand(3)
outband.WriteArray(INTERPMASK)
outband.FlushCache()
outband = outRaster.GetRasterBand(4)
outband.WriteArray(CHIPSIZEX)
outband.FlushCache()
if inps.offset2vx is not None:
ds = gdal.Open(inps.offset2vx)
band = ds.GetRasterBand(1)
offset2vx_1 = band.ReadAsArray()
band = ds.GetRasterBand(2)
offset2vx_2 = band.ReadAsArray()
band=None
ds=None
ds = gdal.Open(inps.offset2vy)
band = ds.GetRasterBand(1)
offset2vy_1 = band.ReadAsArray()
band = ds.GetRasterBand(2)
offset2vy_2 = band.ReadAsArray()
band=None
ds=None
VX = offset2vx_1 * DX + offset2vx_2 * DY
VY = offset2vy_1 * DX + offset2vy_2 * DY
VX = VX.astype(np.float32)
VY = VY.astype(np.float32)
############ write velocity output in Geotiff format
outRaster = driver.Create("velocity.tif", int(xGrid.shape[1]), int(xGrid.shape[0]), 2, gdal.GDT_Float32)
outRaster.SetGeoTransform(tran)
outRaster.SetProjection(proj)
outband = outRaster.GetRasterBand(1)
outband.WriteArray(VX)
outband.FlushCache()
outband = outRaster.GetRasterBand(2)
outband.WriteArray(VY)
outband.FlushCache()
############ prepare for netCDF packaging
if inps.nc_sensor is not None:
vxrefname = str.split(runCmd('fgrep "Velocities:" testGeogrid.txt'))[1]
vyrefname = str.split(runCmd('fgrep "Velocities:" testGeogrid.txt'))[2]
sxname = str.split(runCmd('fgrep "Slopes:" testGeogrid.txt'))[1][:-4]+"s.tif"
syname = str.split(runCmd('fgrep "Slopes:" testGeogrid.txt'))[2][:-4]+"s.tif"
maskname = str.split(runCmd('fgrep "Slopes:" testGeogrid.txt'))[2][:-8]+"masks.tif"
xoff = int(str.split(runCmd('fgrep "Origin index (in DEM) of geogrid:" testGeogrid.txt'))[6])
yoff = int(str.split(runCmd('fgrep "Origin index (in DEM) of geogrid:" testGeogrid.txt'))[7])
xcount = int(str.split(runCmd('fgrep "Dimensions of geogrid:" testGeogrid.txt'))[3])
ycount = int(str.split(runCmd('fgrep "Dimensions of geogrid:" testGeogrid.txt'))[5])
if urlflag == 1:
ds = gdal.Open('/vsicurl/%s' %(vxrefname))
else:
ds = gdal.Open(vxrefname)
band = ds.GetRasterBand(1)
VXref = band.ReadAsArray(xoff, yoff, xcount, ycount)
ds = None
band = None
if urlflag == 1:
ds = gdal.Open('/vsicurl/%s' %(vyrefname))
else:
ds = gdal.Open(vyrefname)
band = ds.GetRasterBand(1)
VYref = band.ReadAsArray(xoff, yoff, xcount, ycount)
ds = None
band = None
if urlflag == 1:
ds = gdal.Open('/vsicurl/%s' %(sxname))
else:
ds = gdal.Open(sxname)
band = ds.GetRasterBand(1)
SX = band.ReadAsArray(xoff, yoff, xcount, ycount)
ds = None
band = None
if urlflag == 1:
ds = gdal.Open('/vsicurl/%s' %(syname))
else:
ds = gdal.Open(syname)
band = ds.GetRasterBand(1)
SY = band.ReadAsArray(xoff, yoff, xcount, ycount)
ds = None
band = None
if urlflag == 1:
ds = gdal.Open('/vsicurl/%s' %(maskname))
else:
ds = gdal.Open(maskname)
band = ds.GetRasterBand(1)
MM = band.ReadAsArray(xoff, yoff, xcount, ycount)
ds = None
band = None
DXref = offset2vy_2 / (offset2vx_1 * offset2vy_2 - offset2vx_2 * offset2vy_1) * VXref - offset2vx_2 / (offset2vx_1 * offset2vy_2 - offset2vx_2 * offset2vy_1) * VYref
DYref = offset2vx_1 / (offset2vx_1 * offset2vy_2 - offset2vx_2 * offset2vy_1) * VYref - offset2vy_1 / (offset2vx_1 * offset2vy_2 - offset2vx_2 * offset2vy_1) * VXref
stable_count = np.sum(SSM & np.logical_not(np.isnan(DX)) & (DX-DXref > -5) & (DX-DXref < 5) & (DY-DYref > -5) & (DY-DYref < 5))
if stable_count == 0:
stable_shift_applied = 0
else:
stable_shift_applied = 1
if stable_shift_applied == 1:
temp = DX.copy() - DXref.copy()
temp[np.logical_not(SSM)] = np.nan
dx_mean_shift = np.median(temp[(temp > -5)&(temp < 5)])
DX = DX - dx_mean_shift
temp = DY.copy() - DYref.copy()
temp[np.logical_not(SSM)] = np.nan
dy_mean_shift = np.median(temp[(temp > -5)&(temp < 5)])
DY = DY - dy_mean_shift
else:
dx_mean_shift = 0.0
dy_mean_shift = 0.0
VX = offset2vx_1 * DX + offset2vx_2 * DY
VY = offset2vy_1 * DX + offset2vy_2 * DY
VX = VX.astype(np.float32)
VY = VY.astype(np.float32)
########################################################################################
############ netCDF packaging for Sentinel and Landsat dataset; can add other sensor format as well
if inps.nc_sensor == "S":
rangePixelSize = float(str.split(runCmd('fgrep "Ground range pixel size:" testGeogrid.txt'))[4])
azimuthPixelSize = float(str.split(runCmd('fgrep "Azimuth pixel size:" testGeogrid.txt'))[3])
dt = float(str.split(runCmd('fgrep "Repeat Time:" testGeogrid.txt'))[2])
epsg = float(str.split(runCmd('fgrep "EPSG:" testGeogrid.txt'))[1])
# print (str(rangePixelSize)+" "+str(azimuthPixelSize))
runCmd('topsinsar_filename.py')
# import scipy.io as sio
conts = sio.loadmat('topsinsar_filename.mat')
master_filename = conts['master_filename'][0]
slave_filename = conts['slave_filename'][0]
master_dt = conts['master_dt'][0]
slave_dt = conts['slave_dt'][0]
master_split = str.split(master_filename,'_')
slave_split = str.split(slave_filename,'_')
import netcdf_output as no
version = '1.0.7'
pair_type = 'radar'
detection_method = 'feature'
coordinates = 'radar'
roi_valid_percentage = int(round(np.sum(CHIPSIZEX!=0)/np.sum(SEARCHLIMITX!=0)*1000.0))/1000
# out_nc_filename = 'Jakobshavn.nc'
PPP = roi_valid_percentage * 100
out_nc_filename = master_filename[0:-4]+'_X_'+slave_filename[0:-4]+'_G0240V02_P{0}{1}{2}.nc'.format(int(PPP/10),int((PPP-int(PPP/10)*10)),int((PPP-int(PPP/10)*10-int((PPP-int(PPP/10)*10)))*10))
out_nc_filename = './' + out_nc_filename
CHIPSIZEY = np.round(CHIPSIZEX * ScaleChipSizeY / 2) * 2
from datetime import date
d0 = date(np.int(master_split[5][0:4]),np.int(master_split[5][4:6]),np.int(master_split[5][6:8]))
d1 = date(np.int(slave_split[5][0:4]),np.int(slave_split[5][4:6]),np.int(slave_split[5][6:8]))
date_dt_base = d1 - d0
date_dt = np.float64(np.abs(date_dt_base.days))
if date_dt_base.days < 0:
date_ct = d1 + (d0 - d1)/2
date_center = date_ct.strftime("%Y%m%d")
else:
date_ct = d0 + (d1 - d0)/2
date_center = date_ct.strftime("%Y%m%d")
IMG_INFO_DICT = {'mission_img1':master_split[0][0],'sensor_img1':'C','satellite_img1':master_split[0][1:3],'acquisition_img1':master_dt,'absolute_orbit_number_img1':master_split[7],'mission_data_take_ID_img1':master_split[8],'product_unique_ID_img1':master_split[9][0:4],'mission_img2':slave_split[0][0],'sensor_img2':'C','satellite_img2':slave_split[0][1:3],'acquisition_img2':slave_dt,'absolute_orbit_number_img2':slave_split[7],'mission_data_take_ID_img2':slave_split[8],'product_unique_ID_img2':slave_split[9][0:4],'date_dt':date_dt,'date_center':date_center,'roi_valid_percentage':roi_valid_percentage,'autoRIFT_software_version':version}
error_vector = np.array([[0.0356, 0.0501, 0.0266, 0.0622, 0.0357, 0.0501],
[0.5194, 1.1638, 0.3319, 1.3701, 0.5191, 1.1628]])
no.netCDF_packaging(VX, VY, DX, DY, INTERPMASK, CHIPSIZEX, CHIPSIZEY, SSM, SX, SY, offset2vx_1, offset2vx_2, offset2vy_1, offset2vy_2, MM, VXref, VYref, rangePixelSize, azimuthPixelSize, dt, epsg, srs, tran, out_nc_filename, pair_type, detection_method, coordinates, IMG_INFO_DICT, stable_count, stable_shift_applied, dx_mean_shift, dy_mean_shift, error_vector)
elif inps.nc_sensor == "L":
XPixelSize = float(str.split(runCmd('fgrep "X-direction pixel size:" testGeogrid.txt'))[3])
YPixelSize = float(str.split(runCmd('fgrep "Y-direction pixel size:" testGeogrid.txt'))[3])
epsg = float(str.split(runCmd('fgrep "EPSG:" testGeogrid.txt'))[1])
master_path = inps.indir_m
slave_path = inps.indir_s
import os
master_filename = os.path.basename(master_path)
slave_filename = os.path.basename(slave_path)
master_split = str.split(master_filename,'_')
slave_split = str.split(slave_filename,'_')
# master_MTL_path = master_path[:-6]+'MTL.txt'
# slave_MTL_path = slave_path[:-6]+'MTL.txt'
#
# master_time = str.split(str.split(runCmd('fgrep "SCENE_CENTER_TIME" '+master_MTL_path))[2][1:-2],':')
# slave_time = str.split(str.split(runCmd('fgrep "SCENE_CENTER_TIME" '+slave_MTL_path))[2][1:-2],':')
master_time = ['00','00','00']
slave_time = ['00','00','00']
from datetime import time as time1
master_time = time1(int(master_time[0]),int(master_time[1]),int(float(master_time[2])))
slave_time = time1(int(slave_time[0]),int(slave_time[1]),int(float(slave_time[2])))
import netcdf_output as no
version = '1.0.7'
pair_type = 'optical'
detection_method = 'feature'
coordinates = 'map'
roi_valid_percentage = int(round(np.sum(CHIPSIZEX!=0)/np.sum(SEARCHLIMITX!=0)*1000.0))/1000
# out_nc_filename = 'Jakobshavn_opt.nc'
PPP = roi_valid_percentage * 100
out_nc_filename = master_filename[0:-8]+'_X_'+slave_filename[0:-8]+'_G0240V02_P{0}{1}{2}.nc'.format(int(PPP/10),int((PPP-int(PPP/10)*10)),int((PPP-int(PPP/10)*10-int((PPP-int(PPP/10)*10)))*10))
out_nc_filename = './' + out_nc_filename
CHIPSIZEY = np.round(CHIPSIZEX * ScaleChipSizeY / 2) * 2
from datetime import date
d0 = date(np.int(master_split[3][0:4]),np.int(master_split[3][4:6]),np.int(master_split[3][6:8]))
d1 = date(np.int(slave_split[3][0:4]),np.int(slave_split[3][4:6]),np.int(slave_split[3][6:8]))
date_dt_base = d1 - d0
date_dt = np.float64(np.abs(date_dt_base.days))
if date_dt_base.days < 0:
date_ct = d1 + (d0 - d1)/2
date_center = date_ct.strftime("%Y%m%d")
else:
date_ct = d0 + (d1 - d0)/2
date_center = date_ct.strftime("%Y%m%d")
master_dt = master_split[3][0:8] + master_time.strftime("T%H:%M:%S")
slave_dt = slave_split[3][0:8] + slave_time.strftime("T%H:%M:%S")
IMG_INFO_DICT = {'mission_img1':master_split[0][0],'sensor_img1':master_split[0][1],'satellite_img1':np.float64(master_split[0][2:4]),'correction_level_img1':master_split[1],'path_img1':np.float64(master_split[2][0:3]),'row_img1':np.float64(master_split[2][3:6]),'acquisition_date_img1':master_dt,'processing_date_img1':master_split[4][0:8],'collection_number_img1':np.float64(master_split[5]),'collection_category_img1':master_split[6],'mission_img2':slave_split[0][0],'sensor_img2':slave_split[0][1],'satellite_img2':np.float64(slave_split[0][2:4]),'correction_level_img2':slave_split[1],'path_img2':np.float64(slave_split[2][0:3]),'row_img2':np.float64(slave_split[2][3:6]),'acquisition_date_img2':slave_dt,'processing_date_img2':slave_split[4][0:8],'collection_number_img2':np.float64(slave_split[5]),'collection_category_img2':slave_split[6],'date_dt':date_dt,'date_center':date_center,'roi_valid_percentage':roi_valid_percentage,'autoRIFT_software_version':version}
error_vector = np.array([57.,57.])
no.netCDF_packaging(VX, VY, DX, DY, INTERPMASK, CHIPSIZEX, CHIPSIZEY, SSM, SX, SY, offset2vx_1, offset2vx_2, offset2vy_1, offset2vy_2, MM, VXref, VYref, XPixelSize, YPixelSize, None, epsg, srs, tran, out_nc_filename, pair_type, detection_method, coordinates, IMG_INFO_DICT, stable_count, stable_shift_applied, dx_mean_shift, dy_mean_shift, error_vector)
elif inps.nc_sensor == "S2":
XPixelSize = float(str.split(runCmd('fgrep "X-direction pixel size:" testGeogrid.txt'))[3])
YPixelSize = float(str.split(runCmd('fgrep "Y-direction pixel size:" testGeogrid.txt'))[3])
epsg = float(str.split(runCmd('fgrep "EPSG:" testGeogrid.txt'))[1])
master_path = inps.indir_m
slave_path = inps.indir_s
master_split = master_path.split('_')
slave_split = slave_path.split('_')
import os
master_filename = master_split[0][-3:]+'_'+master_split[2]+'_'+master_split[4][:3]+'_'+os.path.basename(master_path)
slave_filename = slave_split[0][-3:]+'_'+slave_split[2]+'_'+slave_split[4][:3]+'_'+os.path.basename(slave_path)
master_time = ['00','00','00']
slave_time = ['00','00','00']
from datetime import time as time1
master_time = time1(int(master_time[0]),int(master_time[1]),int(float(master_time[2])))
slave_time = time1(int(slave_time[0]),int(slave_time[1]),int(float(slave_time[2])))
import netcdf_output as no
version = '1.0.7'
pair_type = 'optical'
detection_method = 'feature'
coordinates = 'map'
roi_valid_percentage = int(round(np.sum(CHIPSIZEX!=0)/np.sum(SEARCHLIMITX!=0)*1000.0))/1000
PPP = roi_valid_percentage * 100
out_nc_filename = master_filename[0:-8]+'_X_'+slave_filename[0:-8]+'_G0240V02_P{0}{1}{2}.nc'.format(int(PPP/10),int((PPP-int(PPP/10)*10)),int((PPP-int(PPP/10)*10-int((PPP-int(PPP/10)*10)))*10))
out_nc_filename = './' + out_nc_filename
CHIPSIZEY = np.round(CHIPSIZEX * ScaleChipSizeY / 2) * 2
from datetime import date
d0 = date(np.int(master_split[2][0:4]),np.int(master_split[2][4:6]),np.int(master_split[2][6:8]))
d1 = date(np.int(slave_split[2][0:4]),np.int(slave_split[2][4:6]),np.int(slave_split[2][6:8]))
date_dt_base = d1 - d0
date_dt = np.float64(np.abs(date_dt_base.days))
if date_dt_base.days < 0:
date_ct = d1 + (d0 - d1)/2
date_center = date_ct.strftime("%Y%m%d")
else:
date_ct = d0 + (d1 - d0)/2
date_center = date_ct.strftime("%Y%m%d")
master_dt = master_split[2] + master_time.strftime("T%H:%M:%S")
slave_dt = slave_split[2] + slave_time.strftime("T%H:%M:%S")
IMG_INFO_DICT = {'mission_img1':master_split[0][-3],'satellite_img1':master_split[0][-2:],'correction_level_img1':master_split[4][:3],'acquisition_date_img1':master_dt,'mission_img2':slave_split[0][-3],'satellite_img2':slave_split[0][-2:],'correction_level_img2':slave_split[4][:3],'acquisition_date_img2':slave_dt,'date_dt':date_dt,'date_center':date_center,'roi_valid_percentage':roi_valid_percentage,'autoRIFT_software_version':version}
error_vector = np.array([57.,57.])
no.netCDF_packaging(VX, VY, DX, DY, INTERPMASK, CHIPSIZEX, CHIPSIZEY, SSM, SX, SY, offset2vx_1, offset2vx_2, offset2vy_1, offset2vy_2, MM, VXref, VYref, XPixelSize, YPixelSize, None, epsg, srs, tran, out_nc_filename, pair_type, detection_method, coordinates, IMG_INFO_DICT, stable_count, stable_shift_applied, dx_mean_shift, dy_mean_shift, error_vector)
elif inps.nc_sensor is None:
print('netCDF packaging not performed')
else:
raise Exception('netCDF packaging not supported for the type "{0}"'.format(inps.nc_sensor))
print("Write Outputs Done!!!")
print(time.time()-t1)