def select_pairs_sequential(date_list, num_connection=2, date12_format='YYMMDD-YYMMDD'): """Select Pairs in a Sequential way: For each acquisition, find its num_connection nearest acquisitions in the past time. Inputs: date_list : list of date in YYMMDD/YYYYMMDD format Reference: Fattahi, H., and F. Amelung (2013), DEM Error Correction in InSAR Time Series, IEEE TGRS, 51(7), 4249-4259. """ date8_list = sorted(ptime.yyyymmdd(date_list)) date6_list = ptime.yymmdd(date8_list) date_idx_list = list(range(len(date6_list))) # Get pairs index list date12_idx_list = [] for date_idx in date_idx_list: for i in range(num_connection): if date_idx - i - 1 >= 0: date12_idx_list.append([date_idx - i - 1, date_idx]) date12_idx_list = [sorted(idx) for idx in sorted(date12_idx_list)] # Convert index into date12 date12_list = [ date6_list[idx[0]] + '-' + date6_list[idx[1]] for idx in date12_idx_list ] if date12_format == 'YYYYMMDD_YYYYMMDD': date12_list = ptime.yyyymmdd_date12(date12_list) return date12_list
def plot_coherence_matrix(ax, date12List, cohList, date12List_drop=[], plot_dict={}): """Plot Coherence Matrix of input network if date12List_drop is not empty, plot KEPT pairs in the upper triangle and ALL pairs in the lower triangle. """ # Figure Setting if not 'fontsize' in plot_dict.keys(): plot_dict['fontsize'] = 12 if not 'linewidth' in plot_dict.keys(): plot_dict['linewidth'] = 2 if not 'markercolor' in plot_dict.keys(): plot_dict['markercolor'] = 'orange' if not 'markersize' in plot_dict.keys(): plot_dict['markersize'] = 16 if not 'disp_title' in plot_dict.keys(): plot_dict['disp_title'] = True if not 'cbar_label' in plot_dict.keys(): plot_dict['cbar_label'] = 'Coherence' date12List = ptime.yyyymmdd_date12(date12List) coh_mat = pnet.coherence_matrix(date12List, cohList) if date12List_drop: # Date Convert m_dates = [i.split('_')[0] for i in date12List] s_dates = [i.split('_')[1] for i in date12List] dateList = sorted(list(set(m_dates + s_dates))) # Set dropped pairs' value to nan, in upper triangle only. for date12 in date12List_drop: idx1, idx2 = [dateList.index(i) for i in date12.split('_')] coh_mat[idx1, idx2] = np.nan # Show diagonal value as black, to be distinguished from un-selected interferograms diag_mat = np.diag(np.ones(coh_mat.shape[0])) diag_mat[diag_mat == 0.] = np.nan im = ax.imshow(diag_mat, cmap='gray_r', vmin=0.0, vmax=1.0, interpolation='nearest') im = ax.imshow(coh_mat, cmap='jet', vmin=0.0, vmax=1.0, interpolation='nearest') date_num = coh_mat.shape[0] if date_num < 30: tick_list = list(range(0, date_num, 5)) else: tick_list = list(range(0, date_num, 10)) ax.get_xaxis().set_ticks(tick_list) ax.get_yaxis().set_ticks(tick_list) ax.set_xlabel('Image Number', fontsize=plot_dict['fontsize']) ax.set_ylabel('Image Number', fontsize=plot_dict['fontsize']) if plot_dict['disp_title']: ax.set_title('Coherence Matrix') # Colorbar divider = make_axes_locatable(ax) cax = divider.append_axes("right", "3%", pad="3%") cbar = plt.colorbar(im, cax=cax) cbar.set_label(plot_dict['cbar_label'], fontsize=plot_dict['fontsize']) # Legend if date12List_drop: ax.plot([], [], label='Upper: used ifgrams') ax.plot([], [], label='Lower: all ifgrams') ax.legend(handlelength=0) return ax
def select_pairs_all(date_list, date12_format='YYMMDD-YYMMDD'): """Select All Possible Pairs/Interferograms Input : date_list - list of date in YYMMDD/YYYYMMDD format Output: date12_list - list date12 in YYMMDD-YYMMDD format Reference: Berardino, P., G. Fornaro, R. Lanari, and E. Sansosti (2002), A new algorithm for surface deformation monitoring based on small baseline differential SAR interferograms, IEEE TGRS, 40(11), 2375-2383. """ date8_list = sorted(ptime.yyyymmdd(date_list)) date6_list = ptime.yymmdd(date8_list) date12_list = list(itertools.combinations(date6_list, 2)) date12_list = [date12[0] + '-' + date12[1] for date12 in date12_list] if date12_format == 'YYYYMMDD_YYYYMMDD': date12_list = ptime.yyyymmdd_date12(date12_list) return date12_list
def select_pairs_mst(date_list, pbase_list, date12_format='YYMMDD-YYMMDD'): """Select Pairs using Minimum Spanning Tree technique Connection Cost is calculated using the baseline distance in perp and scaled temporal baseline (Pepe and Lanari, 2006, TGRS) plane. Inputs: date_list : list of date in YYMMDD/YYYYMMDD format pbase_list : list of float, perpendicular spatial baseline References: Pepe, A., and R. Lanari (2006), On the extension of the minimum cost flow algorithm for phase unwrapping of multitemporal differential SAR interferograms, IEEE TGRS, 44(9), 2374-2383. Perissin D., Wang T. (2012), Repeat-pass SAR interferometry with partially coherent targets. IEEE TGRS. 271-280 """ # Get temporal baseline in days date6_list = ptime.yymmdd(date_list) date8_list = ptime.yyyymmdd(date_list) tbase_list = ptime.date_list2tbase(date8_list)[0] # Normalization (Pepe and Lanari, 2006, TGRS) temp2perp_scale = (max(pbase_list) - min(pbase_list)) / (max(tbase_list) - min(tbase_list)) tbase_list = [tbase * temp2perp_scale for tbase in tbase_list] # Get weight matrix ttMat1, ttMat2 = np.meshgrid(np.array(tbase_list), np.array(tbase_list)) ppMat1, ppMat2 = np.meshgrid(np.array(pbase_list), np.array(pbase_list)) ttMat = np.abs(ttMat1 - ttMat2) # temporal distance matrix ppMat = np.abs(ppMat1 - ppMat2) # spatial distance matrix # 2D distance matrix in temp/perp domain weightMat = np.sqrt(np.square(ttMat) + np.square(ppMat)) weightMat = sparse.csr_matrix(weightMat) # compress sparse row matrix # MST path based on weight matrix mstMat = sparse.csgraph.minimum_spanning_tree(weightMat) # Convert MST index matrix into date12 list [s_idx_list, m_idx_list] = [ date_idx_array.tolist() for date_idx_array in sparse.find(mstMat)[0:2] ] date12_list = [] for i in range(len(m_idx_list)): idx = sorted([m_idx_list[i], s_idx_list[i]]) date12 = date6_list[idx[0]] + '-' + date6_list[idx[1]] date12_list.append(date12) if date12_format == 'YYYYMMDD_YYYYMMDD': date12_list = ptime.yyyymmdd_date12(date12_list) return date12_list
def select_pairs_star(date_list, m_date=None, pbase_list=[], date12_format='YYMMDD-YYMMDD'): """Select Star-like network/interferograms/pairs, it's a single master network, similar to PS approach. Usage: m_date : master date, choose it based on the following cretiria: 1) near the center in temporal and spatial baseline 2) prefer winter season than summer season for less temporal decorrelation Reference: Ferretti, A., C. Prati, and F. Rocca (2001), Permanent scatterers in SAR interferometry, IEEE TGRS, 39(1), 8-20. """ date8_list = sorted(ptime.yyyymmdd(date_list)) date6_list = ptime.yymmdd(date8_list) # Select master date if not existed if not m_date: m_date = select_master_date(date8_list, pbase_list) print(('auto select master date: ' + m_date)) # Check input master date m_date8 = ptime.yyyymmdd(m_date) if m_date8 not in date8_list: print('Input master date is not existed in date list!') print(('Input master date: ' + str(m_date8))) print(('Input date list: ' + str(date8_list))) m_date8 = None # Generate star/ps network m_idx = date8_list.index(m_date8) date12_idx_list = [ sorted([m_idx, s_idx]) for s_idx in range(len(date8_list)) if s_idx is not m_idx ] date12_list = [ date6_list[idx[0]] + '-' + date6_list[idx[1]] for idx in date12_idx_list ] if date12_format == 'YYYYMMDD_YYYYMMDD': date12_list = ptime.yyyymmdd_date12(date12_list) return date12_list
def select_pairs_delaunay(date_list, pbase_list, norm=True, date12_format='YYMMDD-YYMMDD'): """Select Pairs using Delaunay Triangulation based on temporal/perpendicular baselines Inputs: date_list : list of date in YYMMDD/YYYYMMDD format pbase_list : list of float, perpendicular spatial baseline norm : normalize temporal baseline to perpendicular baseline Key points 1. Define a ratio between perpendicular and temporal baseline axis units (Pepe and Lanari, 2006, TGRS). 2. Pairs with too large perpendicular / temporal baseline or Doppler centroid difference should be removed after this, using a threshold, to avoid strong decorrelations (Zebker and Villasenor, 1992, TGRS). Reference: Pepe, A., and R. Lanari (2006), On the extension of the minimum cost flow algorithm for phase unwrapping of multitemporal differential SAR interferograms, IEEE TGRS, 44(9), 2374-2383. Zebker, H. A., and J. Villasenor (1992), Decorrelation in interferometric radar echoes, IEEE TGRS, 30(5), 950-959. """ # Get temporal baseline in days date6_list = ptime.yymmdd(date_list) date8_list = ptime.yyyymmdd(date_list) tbase_list = ptime.date_list2tbase(date8_list)[0] # Normalization (Pepe and Lanari, 2006, TGRS) if norm: temp2perp_scale = (max(pbase_list) - min(pbase_list)) / ( max(tbase_list) - min(tbase_list)) tbase_list = [tbase * temp2perp_scale for tbase in tbase_list] # Generate Delaunay Triangulation date12_idx_list = Triangulation(tbase_list, pbase_list).edges.tolist() date12_idx_list = [sorted(idx) for idx in sorted(date12_idx_list)] # Convert index into date12 date12_list = [ date6_list[idx[0]] + '-' + date6_list[idx[1]] for idx in date12_idx_list ] if date12_format == 'YYYYMMDD_YYYYMMDD': date12_list = ptime.yyyymmdd_date12(date12_list) return date12_list
def select_pairs_hierarchical(date_list, pbase_list, temp_perp_list, date12_format='YYMMDD-YYMMDD'): """Select Pairs in a hierarchical way using list of temporal and perpendicular baseline thresholds For each temporal/perpendicular combination, select all possible pairs; and then merge all combination results together for the final output (Zhao, 2015). Inputs: date_list : list of date in YYMMDD/YYYYMMDD format pbase_list : list of float, perpendicular spatial baseline temp_perp_list : list of list of 2 floats, for list of temporal/perp baseline, e.g. [[32.0, 800.0], [48.0, 600.0], [64.0, 200.0]] Examples: pairs = select_pairs_hierarchical(date_list, pbase_list, [[32.0, 800.0], [48.0, 600.0], [64.0, 200.0]]) Reference: Zhao, W., (2015), Small deformation detected from InSAR time-series and their applications in geophysics, Doctoral dissertation, Univ. of Miami, Section 6.3. """ # Get all date12 date12_list_all = select_pairs_all(date_list) # Loop of Threshold print('List of temporal and perpendicular spatial baseline thresholds:') print(temp_perp_list) date12_list = [] for temp_perp in temp_perp_list: tbase_max = temp_perp[0] pbase_max = temp_perp[1] date12_list_tmp = threshold_temporal_baseline(date12_list_all, tbase_max, keep_seasonal=False) date12_list_tmp = threshold_perp_baseline(date12_list_tmp, date_list, pbase_list, pbase_max) date12_list += date12_list_tmp date12_list = sorted(list(set(date12_list))) if date12_format == 'YYYYMMDD_YYYYMMDD': date12_list = ptime.yyyymmdd_date12(date12_list) return date12_list
def plot_coherence_history(ax, date12List, cohList, plot_dict={}): """Plot min/max Coherence of all interferograms for each date""" # Figure Setting if not 'fontsize' in plot_dict.keys(): plot_dict['fontsize'] = 12 if not 'linewidth' in plot_dict.keys(): plot_dict['linewidth'] = 2 if not 'markercolor' in plot_dict.keys(): plot_dict['markercolor'] = 'orange' if not 'markersize' in plot_dict.keys(): plot_dict['markersize'] = 16 if not 'disp_title' in plot_dict.keys(): plot_dict['disp_title'] = True if not 'every_year' in plot_dict.keys(): plot_dict['every_year'] = 1 # Get date list date12List = ptime.yyyymmdd_date12(date12List) m_dates = [date12.split('_')[0] for date12 in date12List] s_dates = [date12.split('_')[1] for date12 in date12List] dateList = sorted(ptime.yyyymmdd(list(set(m_dates + s_dates)))) dates, datevector = ptime.date_list2vector(dateList) bar_width = ut.most_common(np.diff(dates).tolist())*3/4 x_list = [i-bar_width/2 for i in dates] coh_mat = pnet.coherence_matrix(date12List, cohList) ax.bar(x_list, np.nanmax(coh_mat, axis=0), bar_width.days, label='Max Coherence') ax.bar(x_list, np.nanmin(coh_mat, axis=0), bar_width.days, label='Min Coherence') if plot_dict['disp_title']: ax.set_title('Coherence History of All Related Interferograms') ax = auto_adjust_xaxis_date(ax, datevector, plot_dict['fontsize'], every_year=plot_dict['every_year'])[0] ax.set_ylim([0.0, 1.0]) ax.set_xlabel('Time [years]', fontsize=plot_dict['fontsize']) ax.set_ylabel('Coherence', fontsize=plot_dict['fontsize']) ax.legend(loc='lower right') return ax
def simulate_coherence(date12_list, baseline_file='bl_list.txt', sensor_name='Env', inc_angle=22.8, decor_time=200.0, coh_resid=0.2, display=False): """Simulate coherence for a given set of interferograms Inputs: date12_list - list of string in YYMMDD-YYMMDD format, indicating pairs configuration baseline_file - string, path of baseline list text file sensor_name - string, SAR sensor name inc_angle - float, incidence angle decor_time - float / 2D np.array in size of (1, pixel_num) decorrelation rate in days, time for coherence to drop to 1/e of its initial value coh_resid - float / 2D np.array in size of (1, pixel_num) long-term coherence, minimum attainable coherence value display - bool, display result as matrix or not Output: cohs - 2D np.array in size of (ifgram_num, pixel_num) Example: date12_list = pnet.get_date12_list('ifgram_list.txt') cohs = simulate_coherences(date12_list, 'bl_list.txt', sensor_name='Tsx') References: Zebker, H. A., & Villasenor, J. (1992). Decorrelation in interferometric radar echoes. IEEE-TGRS, 30(5), 950-959. Hanssen, R. F. (2001). Radar interferometry: data interpretation and error analysis (Vol. 2). Dordrecht, Netherlands: Kluwer Academic Pub. Morishita, Y., & Hanssen, R. F. (2015). Temporal decorrelation in L-, C-, and X-band satellite radar interferometry for pasture on drained peat soils. IEEE-TGRS, 53(2), 1096-1104. Parizzi, A., Cong, X., & Eineder, M. (2009). First Results from Multifrequency Interferometry. A comparison of different decorrelation time constants at L, C, and X Band. ESA Scientific Publications(SP-677), 1-5. """ date_list, pbase_list, dop_list = read_baseline_file(baseline_file)[0:3] tbase_list = ptime.date_list2tbase(date_list)[0] # Thermal decorrelation (Zebker and Villasenor, 1992, Eq.4) SNR = sensor.signal2noise_ratio(sensor_name) coh_thermal = 1. / (1. + 1. / SNR) pbase_c = critical_perp_baseline(sensor_name, inc_angle) bandwidth_az = sensor.azimuth_bandwidth(sensor_name) date12_list = ptime.yyyymmdd_date12(date12_list) ifgram_num = len(date12_list) if isinstance(decor_time, (int, float)): pixel_num = 1 decor_time = float(decor_time) else: pixel_num = decor_time.shape[1] if decor_time == 0.: decor_time = 0.01 cohs = np.zeros((ifgram_num, pixel_num), np.float32) for i in range(ifgram_num): if display: sys.stdout.write('\rinterferogram = %4d/%4d' % (i, ifgram_num)) sys.stdout.flush() m_date, s_date = date12_list[i].split('_') m_idx = date_list.index(m_date) s_idx = date_list.index(s_date) pbase = pbase_list[s_idx] - pbase_list[m_idx] tbase = tbase_list[s_idx] - tbase_list[m_idx] # Geometric decorrelation (Hanssen, 2001, Eq. 4.4.12) coh_geom = (pbase_c - abs(pbase)) / pbase_c if coh_geom < 0.: coh_geom = 0. # Doppler centroid decorrelation (Hanssen, 2001, Eq. 4.4.13) if not dop_list: coh_dc = 1. else: coh_dc = calculate_doppler_overlap(dop_list[m_idx], dop_list[s_idx], bandwidth_az) if coh_dc < 0.: coh_dc = 0. # Option 1: Temporal decorrelation - exponential delay model (Parizzi et al., 2009; Morishita and Hanssen, 2015) coh_temp = np.multiply( (coh_thermal - coh_resid), np.exp( -1 * abs(tbase) / decor_time)) + coh_resid coh = coh_geom * coh_dc * coh_temp cohs[i, :] = coh #epsilon = 1e-3 #cohs[cohs < epsilon] = epsilon if display: print('') if display: print(('critical perp baseline: %.f m' % pbase_c)) cohs_mat = coherence_matrix(date12_list, cohs) plt.figure() plt.imshow(cohs_mat, vmin=0.0, vmax=1.0, cmap='jet') plt.xlabel('Image number') plt.ylabel('Image number') cbar = plt.colorbar() cbar.set_label('Coherence') plt.title('Coherence matrix') plt.show() return cohs
def get_date12_to_drop(inps): """Get date12 list to dropped Return [] if no ifgram to drop, thus KEEP ALL ifgrams; None if nothing to change, exit without doing anything. """ obj = ifgramStack(inps.file) obj.open() date12ListAll = obj.date12List dateList = obj.dateList print('number of interferograms: {}'.format(len(date12ListAll))) # Get date12_to_drop date12_to_drop = [] # reference file if inps.referenceFile: date12_to_keep = ifgramStack(inps.referenceFile).get_date12_list(dropIfgram=True) print('--------------------------------------------------') print('use reference pairs info from file: {}'.format(inps.referenceFile)) print('number of interferograms in reference: {}'.format(len(date12_to_keep))) tempList = sorted(list(set(date12ListAll) - set(date12_to_keep))) date12_to_drop += tempList print('date12 not in reference file: ({})\n{}'.format(len(tempList), tempList)) # coherence file if inps.coherenceBased: print('--------------------------------------------------') print('use coherence-based network modification') coord = ut.coordinate(obj.metadata, lookup_file=inps.lookupFile) if inps.aoi_geo_box and inps.lookupFile: print('input AOI in (lon0, lat1, lon1, lat0): {}'.format(inps.aoi_geo_box)) inps.aoi_pix_box = coord.bbox_geo2radar(inps.aoi_geo_box) if inps.aoi_pix_box: inps.aoi_pix_box = coord.check_box_within_data_coverage(inps.aoi_pix_box) print('input AOI in (x0,y0,x1,y1): {}'.format(inps.aoi_pix_box)) # Calculate spatial average coherence cohList = ut.spatial_average(inps.file, datasetName='coherence', maskFile=inps.maskFile, box=inps.aoi_pix_box, saveList=True)[0] coh_date12_list = list(np.array(date12ListAll)[np.array(cohList) >= inps.minCoherence]) # MST network if inps.keepMinSpanTree: print('Get minimum spanning tree (MST) of interferograms with inverse of coherence.') msg = ('Drop ifgrams with ' '1) average coherence < {} AND ' '2) not in MST network: '.format(inps.minCoherence)) mst_date12_list = pnet.threshold_coherence_based_mst(date12ListAll, cohList) mst_date12_list = ptime.yyyymmdd_date12(mst_date12_list) else: msg = 'Drop ifgrams with average coherence < {}: '.format(inps.minCoherence) mst_date12_list = [] tempList = sorted(list(set(date12ListAll) - set(coh_date12_list + mst_date12_list))) date12_to_drop += tempList msg += '({})'.format(len(tempList)) if len(tempList) <= 200: msg += '\n{}'.format(tempList) print(msg) # temp baseline threshold if inps.tempBaseMax: tempIndex = np.abs(obj.tbaseIfgram) > inps.tempBaseMax tempList = list(np.array(date12ListAll)[tempIndex]) date12_to_drop += tempList print('--------------------------------------------------') print('Drop ifgrams with temporal baseline > {} days: ({})\n{}'.format( inps.tempBaseMax, len(tempList), tempList)) # perp baseline threshold if inps.perpBaseMax: tempIndex = np.abs(obj.pbaseIfgram) > inps.perpBaseMax tempList = list(np.array(date12ListAll)[tempIndex]) date12_to_drop += tempList print('--------------------------------------------------') print('Drop ifgrams with perp baseline > {} meters: ({})\n{}'.format( inps.perpBaseMax, len(tempList), tempList)) # connection number threshold if inps.connNumMax: seq_date12_list = pnet.select_pairs_sequential(dateList, inps.connNumMax) seq_date12_list = ptime.yyyymmdd_date12(seq_date12_list) tempList = [i for i in date12ListAll if i not in seq_date12_list] date12_to_drop += tempList print('--------------------------------------------------') msg = 'Drop ifgrams with temporal baseline beyond {} neighbors: ({})'.format( inps.connNumMax, len(tempList)) if len(tempList) <= 200: msg += '\n{}'.format(tempList) print(msg) # excludeIfgIndex if inps.excludeIfgIndex: tempList = [date12ListAll[i] for i in inps.excludeIfgIndex] date12_to_drop += tempList print('--------------------------------------------------') print('Drop ifgrams with the following index number: {}'.format(len(tempList))) for i in range(len(tempList)): print('{} : {}'.format(i, tempList[i])) #len(tempList), zip(inps.excludeIfgIndex, tempList))) # excludeDate if inps.excludeDate: tempList = [i for i in date12ListAll if any(j in inps.excludeDate for j in i.split('_'))] date12_to_drop += tempList print('-'*50+'\nDrop ifgrams including the following dates: ({})\n{}'.format( len(tempList), inps.excludeDate)) print('-'*30+'\n{}'.format(tempList)) # startDate if inps.startDate: minDate = int(inps.startDate) tempList = [i for i in date12ListAll if any(int(j) < minDate for j in i.split('_'))] date12_to_drop += tempList print('--------------------------------------------------') print('Drop ifgrams with date earlier than: {} ({})\n{}'.format( inps.startDate, len(tempList), tempList)) # endDate if inps.endDate: maxDate = int(inps.endDate) tempList = [i for i in date12ListAll if any(int(j) > maxDate for j in i.split('_'))] date12_to_drop += tempList print('--------------------------------------------------') print('Drop ifgrams with date later than: {} ({})\n{}'.format( inps.endDate, len(tempList), tempList)) # Manually drop pairs if inps.manual: tempList = manual_select_pairs_to_remove(inps.file) if tempList is None: return None tempList = [i for i in tempList if i in date12ListAll] print('date12 selected to remove: ({})\n{}'.format(len(tempList), tempList)) date12_to_drop += tempList # drop duplicate date12 and sort in order date12_to_drop = sorted(list(set(date12_to_drop))) date12_to_keep = sorted(list(set(date12ListAll) - set(date12_to_drop))) print('--------------------------------------------------') print('number of interferograms to remove: {}'.format(len(date12_to_drop))) print('number of interferograms to keep : {}'.format(len(date12_to_keep))) date_to_keep = [d for date12 in date12_to_keep for d in date12.split('_')] date_to_keep = sorted(list(set(date_to_keep))) date_to_drop = sorted(list(set(dateList) - set(date_to_keep))) if len(date_to_drop) > 0: print('number of acquisitions to remove: {}\n{}'.format(len(date_to_drop), date_to_drop)) date12ListKept = obj.get_date12_list(dropIfgram=True) date12ListDropped = sorted(list(set(date12ListAll) - set(date12ListKept))) if date12_to_drop == date12ListDropped: print('Calculated date12 to drop is the same as exsiting marked input file, skip updating file.') date12_to_drop = None elif date12_to_drop == date12ListAll: raise Exception('Zero interferogram left! Please adjust your setting and try again.') return date12_to_drop
def plot_network(ax, date12List, dateList, pbaseList, plot_dict={}, date12List_drop=[], print_msg=True): """Plot Temporal-Perp baseline Network Inputs ax : matplotlib axes object date12List : list of string for date12 in YYYYMMDD_YYYYMMDD format dateList : list of string, for date in YYYYMMDD format pbaseList : list of float, perp baseline, len=number of acquisition plot_dict : dictionary with the following items: fontsize linewidth markercolor markersize cohList : list of float, coherence value of each interferogram, len = number of ifgrams disp_min/max : float, min/max range of the color display based on cohList colormap : string, colormap name coh_thres : float, coherence of where to cut the colormap for display disp_title : bool, show figure title or not, default: True disp_drop: bool, show dropped interferograms or not, default: True Output ax : matplotlib axes object """ # Figure Setting if not 'fontsize' in plot_dict.keys(): plot_dict['fontsize'] = 12 if not 'linewidth' in plot_dict.keys(): plot_dict['linewidth'] = 2 if not 'markercolor' in plot_dict.keys(): plot_dict['markercolor'] = 'orange' if not 'markersize' in plot_dict.keys(): plot_dict['markersize'] = 16 # For colorful display of coherence if not 'cohList' in plot_dict.keys(): plot_dict['cohList'] = None if not 'cbar_label' in plot_dict.keys(): plot_dict['cbar_label'] = 'Average Spatial Coherence' if not 'disp_min' in plot_dict.keys(): plot_dict['disp_min'] = 0.2 if not 'disp_max' in plot_dict.keys(): plot_dict['disp_max'] = 1.0 if not 'colormap' in plot_dict.keys(): plot_dict['colormap'] = 'RdBu' if not 'disp_title' in plot_dict.keys(): plot_dict['disp_title'] = True if not 'coh_thres' in plot_dict.keys(): plot_dict['coh_thres'] = None if not 'disp_drop' in plot_dict.keys(): plot_dict['disp_drop'] = True if not 'every_year' in plot_dict.keys(): plot_dict['every_year'] = 1 cohList = plot_dict['cohList'] disp_min = plot_dict['disp_min'] disp_max = plot_dict['disp_max'] coh_thres = plot_dict['coh_thres'] transparency = 0.7 # Date Convert dateList = ptime.yyyymmdd(sorted(dateList)) dates, datevector = ptime.date_list2vector(dateList) tbaseList = ptime.date_list2tbase(dateList)[0] ## maxBperp and maxBtemp date12List = ptime.yyyymmdd_date12(date12List) ifgram_num = len(date12List) pbase12 = np.zeros(ifgram_num) tbase12 = np.zeros(ifgram_num) for i in range(ifgram_num): m_date, s_date = date12List[i].split('_') m_idx = dateList.index(m_date) s_idx = dateList.index(s_date) pbase12[i] = pbaseList[s_idx] - pbaseList[m_idx] tbase12[i] = tbaseList[s_idx] - tbaseList[m_idx] if print_msg: print('max perpendicular baseline: {:.2f} m'.format(np.max(np.abs(pbase12)))) print('max temporal baseline: {} days'.format(np.max(tbase12))) ## Keep/Drop - date12 date12List_keep = sorted(list(set(date12List) - set(date12List_drop))) idx_date12_keep = [date12List.index(i) for i in date12List_keep] idx_date12_drop = [date12List.index(i) for i in date12List_drop] if not date12List_drop: plot_dict['disp_drop'] = False ## Keep/Drop - date m_dates = [i.split('_')[0] for i in date12List_keep] s_dates = [i.split('_')[1] for i in date12List_keep] dateList_keep = ptime.yyyymmdd(sorted(list(set(m_dates + s_dates)))) dateList_drop = sorted(list(set(dateList) - set(dateList_keep))) idx_date_keep = [dateList.index(i) for i in dateList_keep] idx_date_drop = [dateList.index(i) for i in dateList_drop] # Ploting # ax=fig.add_subplot(111) # Colorbar when conherence is colored if cohList is not None: data_min = min(cohList) data_max = max(cohList) # Normalize normalization = False if normalization: cohList = [(coh-data_min) / (data_min-data_min) for coh in cohList] disp_min = data_min disp_max = data_max if print_msg: print('showing coherence') print(('colormap: '+plot_dict['colormap'])) print(('display range: '+str([disp_min, disp_max]))) print(('data range: '+str([data_min, data_max]))) splitColormap = True if splitColormap: # Use lower/upper part of colormap to emphasis dropped interferograms if not coh_thres: # Find proper cut percentage so that all keep pairs are blue and drop pairs are red cohList_keep = [cohList[i] for i in idx_date12_keep] cohList_drop = [cohList[i] for i in idx_date12_drop] if cohList_drop: coh_thres = max(cohList_drop) else: coh_thres = min(cohList_keep) if coh_thres < disp_min: disp_min = 0.0 if print_msg: print('data range exceed orginal display range, set new display range to: [0.0, %f]' % (disp_max)) c1_num = np.ceil(200.0 * (coh_thres - disp_min) / (disp_max - disp_min)).astype('int') coh_thres = c1_num / 200.0 * (disp_max-disp_min) + disp_min cmap = plt.get_cmap(plot_dict['colormap']) colors1 = cmap(np.linspace(0.0, 0.3, c1_num)) colors2 = cmap(np.linspace(0.6, 1.0, 200 - c1_num)) cmap = LinearSegmentedColormap.from_list('truncate_RdBu', np.vstack((colors1, colors2))) if print_msg: print(('color jump at '+str(coh_thres))) else: cmap = plt.get_cmap(plot_dict['colormap']) divider = make_axes_locatable(ax) cax = divider.append_axes("right", "3%", pad="3%") norm = mpl.colors.Normalize(vmin=disp_min, vmax=disp_max) cbar = mpl.colorbar.ColorbarBase(cax, cmap=cmap, norm=norm) cbar.set_label(plot_dict['cbar_label'], fontsize=plot_dict['fontsize']) # plot low coherent ifgram first and high coherence ifgram later cohList_keep = [cohList[date12List.index(i)] for i in date12List_keep] date12List_keep = [x for _, x in sorted(zip(cohList_keep, date12List_keep))] # Dot - SAR Acquisition if idx_date_keep: x_list = [dates[i] for i in idx_date_keep] y_list = [pbaseList[i] for i in idx_date_keep] ax.plot(x_list, y_list, 'ko', alpha=0.7, ms=plot_dict['markersize'], mfc=plot_dict['markercolor']) if idx_date_drop: x_list = [dates[i] for i in idx_date_drop] y_list = [pbaseList[i] for i in idx_date_drop] ax.plot(x_list, y_list, 'ko', alpha=0.7, ms=plot_dict['markersize'], mfc='gray') ## Line - Pair/Interferogram # interferograms dropped if plot_dict['disp_drop']: for date12 in date12List_drop: date1, date2 = date12.split('_') idx1 = dateList.index(date1) idx2 = dateList.index(date2) x = np.array([dates[idx1], dates[idx2]]) y = np.array([pbaseList[idx1], pbaseList[idx2]]) if cohList is not None: coh = cohList[date12List.index(date12)] coh_idx = (coh - disp_min) / (disp_max - disp_min) ax.plot(x, y, '--', lw=plot_dict['linewidth'], alpha=transparency, c=cmap(coh_idx)) else: ax.plot(x, y, '--', lw=plot_dict['linewidth'], alpha=transparency, c='k') # interferograms kept for date12 in date12List_keep: date1, date2 = date12.split('_') idx1 = dateList.index(date1) idx2 = dateList.index(date2) x = np.array([dates[idx1], dates[idx2]]) y = np.array([pbaseList[idx1], pbaseList[idx2]]) if cohList is not None: coh = cohList[date12List.index(date12)] coh_idx = (coh - disp_min) / (disp_max - disp_min) ax.plot(x, y, '-', lw=plot_dict['linewidth'], alpha=transparency, c=cmap(coh_idx)) else: ax.plot(x, y, '-', lw=plot_dict['linewidth'], alpha=transparency, c='k') if plot_dict['disp_title']: ax.set_title('Interferogram Network', fontsize=plot_dict['fontsize']) # axis format ax = auto_adjust_xaxis_date(ax, datevector, plot_dict['fontsize'], every_year=plot_dict['every_year'])[0] ax = auto_adjust_yaxis(ax, pbaseList, plot_dict['fontsize']) ax.set_xlabel('Time [years]', fontsize=plot_dict['fontsize']) ax.set_ylabel('Perp Baseline [m]', fontsize=plot_dict['fontsize']) # Legend if plot_dict['disp_drop']: solid_line = mlines.Line2D([], [], color='k', ls='solid', label='Interferograms') dash_line = mlines.Line2D([], [], color='k', ls='dashed', label='Interferograms dropped') ax.legend(handles=[solid_line, dash_line]) return ax