|
| 1 | +#!/usr/bin/python3 |
| 2 | + |
| 3 | +# Reader for Lutan-1 SLC data |
| 4 | +# Used Sentinel1.py and ALOS.py as templates |
| 5 | +# Author: Bryan Marfito, EOS-RS |
| 6 | + |
| 7 | + |
| 8 | +import os |
| 9 | +import glob |
| 10 | +import numpy as np |
| 11 | +import xml.etree.ElementTree as ET |
| 12 | +import datetime |
| 13 | +import isce |
| 14 | +import isceobj |
| 15 | +from isceobj.Planet.Planet import Planet |
| 16 | +from iscesys.Component.Component import Component |
| 17 | +from isceobj.Sensor.Sensor import Sensor |
| 18 | +from isceobj.Scene.Frame import Frame |
| 19 | +from isceobj.Orbit.Orbit import StateVector, Orbit |
| 20 | +from isceobj.Planet.AstronomicalHandbook import Const |
| 21 | +from iscesys.DateTimeUtil.DateTimeUtil import DateTimeUtil as DTUtil |
| 22 | +from isceobj.Orbit.OrbitExtender import OrbitExtender |
| 23 | +from osgeo import gdal |
| 24 | +import warnings |
| 25 | + |
| 26 | +lookMap = { 'RIGHT' : -1, |
| 27 | + 'LEFT' : 1} |
| 28 | + |
| 29 | +# Antenna dimensions 9.8 x 3.4 m |
| 30 | +antennaLength = 9.8 |
| 31 | + |
| 32 | +XML = Component.Parameter('xml', |
| 33 | + public_name = 'xml', |
| 34 | + default = None, |
| 35 | + type = str, |
| 36 | + doc = 'Input XML file') |
| 37 | + |
| 38 | + |
| 39 | +TIFF = Component.Parameter('tiff', |
| 40 | + public_name ='tiff', |
| 41 | + default = None, |
| 42 | + type=str, |
| 43 | + doc = 'Input image file') |
| 44 | + |
| 45 | +ORBIT_FILE = Component.Parameter('orbitFile', |
| 46 | + public_name ='orbitFile', |
| 47 | + default = None, |
| 48 | + type=str, |
| 49 | + doc = 'Orbit file') |
| 50 | + |
| 51 | + |
| 52 | +class Lutan1(Sensor): |
| 53 | + |
| 54 | + "Class for Lutan-1 SLC data" |
| 55 | + |
| 56 | + family = 'l1sm' |
| 57 | + logging_name = 'isce.sensor.Lutan1' |
| 58 | + |
| 59 | + parameter_list = (TIFF, ORBIT_FILE) + Sensor.parameter_list |
| 60 | + |
| 61 | + def __init__(self, name = ''): |
| 62 | + super(Lutan1,self).__init__(self.__class__.family, name=name) |
| 63 | + self.frame = Frame() |
| 64 | + self.frame.configure() |
| 65 | + self._xml_root = None |
| 66 | + self.doppler_coeff = None |
| 67 | + |
| 68 | + def parse(self): |
| 69 | + xmlFileName = self.tiff[:-4] + "meta.xml" |
| 70 | + self.xml = xmlFileName |
| 71 | + |
| 72 | + with open(self.xml, 'r') as fid: |
| 73 | + xmlstr = fid.read() |
| 74 | + |
| 75 | + self._xml_root = ET.fromstring(xmlstr) |
| 76 | + self.populateMetadata() |
| 77 | + fid.close() |
| 78 | + |
| 79 | + if self.orbitFile: |
| 80 | + # Check if orbit file exists or not |
| 81 | + if os.path.isfile(self.orbitFile) == True: |
| 82 | + orb = self.extractOrbit() |
| 83 | + self.frame.orbit.setOrbitSource(os.path.basename(self.orbitFile)) |
| 84 | + else: |
| 85 | + pass |
| 86 | + else: |
| 87 | + warnings.warn("WARNING! No orbit file found. Orbit information from the annotation file is used for processing.") |
| 88 | + orb = self.extractOrbitFromAnnotation() |
| 89 | + self.frame.orbit.setOrbitSource(os.path.basename(self.xml)) |
| 90 | + self.frame.orbit.setOrbitSource('Annotation') |
| 91 | + |
| 92 | + for sv in orb: |
| 93 | + self.frame.orbit.addStateVector(sv) |
| 94 | + |
| 95 | + def convertToDateTime(self,string): |
| 96 | + dt = datetime.datetime.strptime(string,"%Y-%m-%dT%H:%M:%S.%f") |
| 97 | + return dt |
| 98 | + |
| 99 | + |
| 100 | + def grab_from_xml(self, path): |
| 101 | + try: |
| 102 | + res = self._xml_root.find(path).text |
| 103 | + except: |
| 104 | + raise Exception('Tag= %s not found'%(path)) |
| 105 | + |
| 106 | + if res is None: |
| 107 | + raise Exception('Tag = %s not found'%(path)) |
| 108 | + |
| 109 | + return res |
| 110 | + |
| 111 | + |
| 112 | + def populateMetadata(self): |
| 113 | + mission = self.grab_from_xml('generalHeader/mission') |
| 114 | + polarization = self.grab_from_xml('productInfo/acquisitionInfo/polarisationMode') |
| 115 | + frequency = float(self.grab_from_xml('instrument/radarParameters/centerFrequency')) |
| 116 | + passDirection = self.grab_from_xml('productInfo/missionInfo/orbitDirection') |
| 117 | + rangePixelSize = float(self.grab_from_xml('productInfo/imageDataInfo/imageRaster/columnSpacing')) |
| 118 | + azimuthPixelSize = float(self.grab_from_xml('productInfo/imageDataInfo/imageRaster/rowSpacing')) |
| 119 | + rangeSamplingRate = Const.c/(2.0*rangePixelSize) |
| 120 | + |
| 121 | + prf = float(self.grab_from_xml('instrument/settings/settingRecord/PRF')) |
| 122 | + lines = int(self.grab_from_xml('productInfo/imageDataInfo/imageRaster/numberOfRows')) |
| 123 | + samples = int(self.grab_from_xml('productInfo/imageDataInfo/imageRaster/numberOfColumns')) |
| 124 | + |
| 125 | + startingRange = float(self.grab_from_xml('productInfo/sceneInfo/rangeTime/firstPixel'))*Const.c/2.0 |
| 126 | + #slantRange = float(self.grab_from_xml('productSpecific/complexImageInfo/')) |
| 127 | + incidenceAngle = float(self.grab_from_xml('productInfo/sceneInfo/sceneCenterCoord/incidenceAngle')) |
| 128 | + dataStartTime = self.convertToDateTime(self.grab_from_xml('productInfo/sceneInfo/start/timeUTC')) |
| 129 | + dataStopTime = self.convertToDateTime(self.grab_from_xml('productInfo/sceneInfo/stop/timeUTC')) |
| 130 | + pulseLength = float(self.grab_from_xml('processing/processingParameter/rangeCompression/chirps/referenceChirp/pulseLength')) |
| 131 | + pulseBandwidth = float(self.grab_from_xml('processing/processingParameter/rangeCompression/chirps/referenceChirp/pulseBandwidth')) |
| 132 | + chirpSlope = pulseBandwidth/pulseLength |
| 133 | + |
| 134 | + if self.grab_from_xml('processing/processingParameter/rangeCompression/chirps/referenceChirp/chirpSlope') == "DOWN": |
| 135 | + chirpSlope = -1.0 * chirpSlope |
| 136 | + else: |
| 137 | + pass |
| 138 | + |
| 139 | + # Check for satellite's look direction |
| 140 | + if self.grab_from_xml('productInfo/acquisitionInfo/lookDirection') == "LEFT": |
| 141 | + lookSide = lookMap['LEFT'] |
| 142 | + print("Look direction: LEFT") |
| 143 | + else: |
| 144 | + lookSide = lookMap['RIGHT'] |
| 145 | + print("Look direction: RIGHT") |
| 146 | + |
| 147 | + processingFacility = self.grab_from_xml('productInfo/generationInfo/level1ProcessingFacility') |
| 148 | + |
| 149 | + # Platform parameters |
| 150 | + platform = self.frame.getInstrument().getPlatform() |
| 151 | + platform.setPlanet(Planet(pname='Earth')) |
| 152 | + platform.setMission(mission) |
| 153 | + platform.setPointingDirection(lookSide) |
| 154 | + platform.setAntennaLength(antennaLength) |
| 155 | + |
| 156 | + # Instrument parameters |
| 157 | + instrument = self.frame.getInstrument() |
| 158 | + instrument.setRadarFrequency(frequency) |
| 159 | + instrument.setPulseRepetitionFrequency(prf) |
| 160 | + instrument.setPulseLength(pulseLength) |
| 161 | + instrument.setChirpSlope(chirpSlope) |
| 162 | + instrument.setIncidenceAngle(incidenceAngle) |
| 163 | + instrument.setRangePixelSize(rangePixelSize) |
| 164 | + instrument.setRangeSamplingRate(rangeSamplingRate) |
| 165 | + instrument.setPulseLength(pulseLength) |
| 166 | + |
| 167 | + # Frame parameters |
| 168 | + self.frame.setSensingStart(dataStartTime) |
| 169 | + self.frame.setSensingStop(dataStopTime) |
| 170 | + self.frame.setProcessingFacility(processingFacility) |
| 171 | + |
| 172 | + # Two-way travel time |
| 173 | + diffTime = DTUtil.timeDeltaToSeconds(dataStopTime - dataStartTime) / 2.0 |
| 174 | + sensingMid = dataStartTime + datetime.timedelta(microseconds=int(diffTime*1e6)) |
| 175 | + self.frame.setSensingMid(sensingMid) |
| 176 | + self.frame.setPassDirection(passDirection) |
| 177 | + self.frame.setPolarization(polarization) |
| 178 | + self.frame.setStartingRange(startingRange) |
| 179 | + self.frame.setFarRange(startingRange + (samples - 1) * rangePixelSize) |
| 180 | + self.frame.setNumberOfLines(lines) |
| 181 | + self.frame.setNumberOfSamples(samples) |
| 182 | + |
| 183 | + return |
| 184 | + |
| 185 | + |
| 186 | + def extractOrbit(self): |
| 187 | + |
| 188 | + ''' |
| 189 | + Extract orbit information from the orbit file |
| 190 | + ''' |
| 191 | + |
| 192 | + try: |
| 193 | + fp = open(self.orbitFile, 'r') |
| 194 | + except IOError as strerr: |
| 195 | + print("IOError: %s" % strerr) |
| 196 | + |
| 197 | + _xml_root = ET.ElementTree(file=fp).getroot() |
| 198 | + node = _xml_root.find('Data_Block/List_of_OSVs') |
| 199 | + |
| 200 | + orb = Orbit() |
| 201 | + orb.configure() |
| 202 | + |
| 203 | + # I based the margin on the data that I have. |
| 204 | + # Lutan-1 position and velocity sampling frequency is 1 Hz |
| 205 | + margin = datetime.timedelta(seconds=1.0) |
| 206 | + tstart = self.frame.getSensingStart() - margin |
| 207 | + tend = self.frame.getSensingStop() + margin |
| 208 | + |
| 209 | + for child in node: |
| 210 | + timestamp = self.convertToDateTime(child.find('UTC').text) |
| 211 | + if (timestamp >= tstart) and (timestamp <= tend): |
| 212 | + pos = [] |
| 213 | + vel = [] |
| 214 | + for tag in ['VX', 'VY', 'VZ']: |
| 215 | + vel.append(float(child.find(tag).text)) |
| 216 | + |
| 217 | + for tag in ['X', 'Y', 'Z']: |
| 218 | + pos.append(float(child.find(tag).text)) |
| 219 | + |
| 220 | + vec = StateVector() |
| 221 | + vec.setTime(timestamp) |
| 222 | + vec.setPosition(pos) |
| 223 | + vec.setVelocity(vel) |
| 224 | + orb.addStateVector(vec) |
| 225 | + |
| 226 | + fp.close() |
| 227 | + |
| 228 | + return orb |
| 229 | + |
| 230 | + def extractOrbitFromAnnotation(self): |
| 231 | + |
| 232 | + ''' |
| 233 | + Extract orbit information from xml annotation |
| 234 | + WARNING! Only use this method if orbit file is not available |
| 235 | + ''' |
| 236 | + |
| 237 | + try: |
| 238 | + fp = open(self.xml, 'r') |
| 239 | + except IOError as strerr: |
| 240 | + print("IOError: %s" % strerr) |
| 241 | + |
| 242 | + _xml_root = ET.ElementTree(file=fp).getroot() |
| 243 | + node = _xml_root.find('platform/orbit') |
| 244 | + countNode = len(list(_xml_root.find('platform/orbit'))) |
| 245 | + |
| 246 | + frameOrbit = Orbit() |
| 247 | + frameOrbit.setOrbitSource('Header') |
| 248 | + margin = datetime.timedelta(seconds=1.0) |
| 249 | + tstart = self.frame.getSensingStart() - margin |
| 250 | + tend = self.frame.getSensingStop() + margin |
| 251 | + |
| 252 | + for k in range(1,countNode): |
| 253 | + timestamp = self.convertToDateTime(node.find('stateVec[{}]/timeUTC'.format(k)).text) |
| 254 | + if (timestamp >= tstart) and (timestamp <= tend): |
| 255 | + pos = [float(node.find('stateVec[{}]/posX'.format(k)).text), float(node.find('stateVec[{}]/posY'.format(k)).text), float(node.find('stateVec[{}]/posZ'.format(k)).text)] |
| 256 | + vel = [float(node.find('stateVec[{}]/velX'.format(k)).text), float(node.find('stateVec[{}]/velY'.format(k)).text), float(node.find('stateVec[{}]/velZ'.format(k)).text)] |
| 257 | + |
| 258 | + vec = StateVector() |
| 259 | + vec.setTime(timestamp) |
| 260 | + vec.setPosition(pos) |
| 261 | + vec.setVelocity(vel) |
| 262 | + frameOrbit.addStateVector(vec) |
| 263 | + |
| 264 | + fp.close() |
| 265 | + return frameOrbit |
| 266 | + |
| 267 | + def extractImage(self): |
| 268 | + self.parse() |
| 269 | + width = self.frame.getNumberOfSamples() |
| 270 | + lgth = self.frame.getNumberOfLines() |
| 271 | + src = gdal.Open(self.tiff.strip(), gdal.GA_ReadOnly) |
| 272 | + |
| 273 | + # Band 1 as real and band 2 as imaginary numbers |
| 274 | + # Confirmed by Zhang Yunjun |
| 275 | + band1 = src.GetRasterBand(1) |
| 276 | + band2 = src.GetRasterBand(2) |
| 277 | + cJ = np.complex64(1.0j) |
| 278 | + |
| 279 | + fid = open(self.output, 'wb') |
| 280 | + for ii in range(lgth): |
| 281 | + # Combine the real and imaginary to make |
| 282 | + # them in to complex numbers |
| 283 | + real = band1.ReadAsArray(0,ii,width,1) |
| 284 | + imag = band2.ReadAsArray(0,ii,width,1) |
| 285 | + # Data becomes np.complex128 after combining them |
| 286 | + data = real + (cJ * imag) |
| 287 | + data.tofile(fid) |
| 288 | + |
| 289 | + fid.close() |
| 290 | + real = None |
| 291 | + imag = None |
| 292 | + src = None |
| 293 | + band1 = None |
| 294 | + band2 = None |
| 295 | + |
| 296 | + #### |
| 297 | + slcImage = isceobj.createSlcImage() |
| 298 | + slcImage.setByteOrder('l') |
| 299 | + slcImage.setFilename(self.output) |
| 300 | + slcImage.setAccessMode('read') |
| 301 | + slcImage.setWidth(self.frame.getNumberOfSamples()) |
| 302 | + slcImage.setLength(self.frame.getNumberOfLines()) |
| 303 | + slcImage.setXmin(0) |
| 304 | + slcImage.setXmax(self.frame.getNumberOfSamples()) |
| 305 | + self.frame.setImage(slcImage) |
| 306 | + |
| 307 | + def extractDoppler(self): |
| 308 | + ''' |
| 309 | + Set doppler values to zero since the metadata doppler values are unreliable. |
| 310 | + Also, the SLC images are zero doppler. |
| 311 | + ''' |
| 312 | + dop = [0., 0., 0.] |
| 313 | + |
| 314 | + ####For insarApp |
| 315 | + quadratic = {} |
| 316 | + quadratic['a'] = dop[0] / self.frame.getInstrument().getPulseRepetitionFrequency() |
| 317 | + quadratic['b'] = 0. |
| 318 | + quadratic['c'] = 0. |
| 319 | + |
| 320 | + print("Average doppler: ", dop) |
| 321 | + self.frame._dopplerVsPixel = dop |
| 322 | + |
| 323 | + return quadratic |
0 commit comments