def detect_smile_effect(sensor_dict, atm_lut_file): """ Detect smile effect. Parameters ---------- sensor_dict: dict Sensor configurations. atm_lut_file: str Raw atmosphere lookup table filename. """ if os.path.exists(sensor_dict['smile_effect_at_atm_features_file'] ) and os.path.exists(sensor_dict['smile_effect_file']): logger.info( 'Write smile effect at atmosphere aborption features to %s.' % sensor_dict['smile_effect_at_atm_features_file']) logger.info('Write smile effect to %s.' % sensor_dict['smile_effect_file']) return from ENVI import empty_envi_header, read_envi_header, write_envi_header from Spectra import get_closest_wave from scipy import optimize, interpolate import json # Read averaged radiance header = read_envi_header(sensor_dict['avg_rdn_file'] + '.hdr') sensor_rdn = np.memmap(sensor_dict['avg_rdn_file'], mode='r', dtype='float32', shape=(header['lines'], header['samples'])) # shape=(bands, samples) samples = header['samples'] sensor_wave = np.array([float(v) for v in header['waves'].split(',')]) sensor_fwhm = np.array([float(v) for v in header['fwhms'].split(',')]) tmp_vza = header['VZA'].split(',') tmp_raa = header['RAA'].split(',') vza = [] raa = [] for i in range(samples): try: vza.append(float(tmp_vza[i])) raa.append(float(tmp_raa[i])) except: vza.append(np.nan) raa.append(np.nan) logger.info('No spectrum for column: %s.' % (i + 1)) del tmp_vza, tmp_raa del header vza = np.array(vza) raa = np.array(raa) if np.all(np.isnan(vza)) or np.all(np.isnan(raa)): raise IOError( "Cannot detect smile effects since all columns do not have spectra." ) nonnan_index = np.where((~np.isnan(vza)) & (~np.isnan(vza)))[0] vza = np.interp(np.arange(samples), nonnan_index, vza[nonnan_index]) raa = np.interp(np.arange(samples), nonnan_index, raa[nonnan_index]) # Assign visibility vis = [40] * samples # Estimate water vapor column logger.info('Estimate WVC from averaged radiance spectra.') wvc_models = json.load(open(sensor_dict['wvc_model_file'], 'r')) wvc = np.zeros(samples) for model in wvc_models.values(): ratio = sensor_rdn[model['band'][1], :] / ( sensor_rdn[model['band'][0], :] * model['weight'][0] + sensor_rdn[model['band'][2], :] * model['weight'][1]) wvc += np.interp(ratio, model['ratio'], model['wvc']) del ratio wvc /= len(wvc_models) wvc_isnan = np.isnan(wvc) # Look for NaN values if np.any(wvc_isnan): logger.info('WVC could not be estimated for some columns. ' 'Missing values will be interpolated.') interpolate_values(wvc, wvc_isnan) # Replace NaNs with interpolated values logger.info('WVC [mm] statistics: min=%.2f, max=%.2f, avg=%.2f, sd=%.2f.' % (wvc.min(), wvc.max(), wvc.mean(), wvc.std())) del wvc_models, model # Interpolate atmospheric lookup table logger.info('Interpolate atmospheric look-up table.') atm_lut_wave, atm_lut_rdn = interp_atm_lut(atm_lut_file, wvc, vis, vza, raa) # shape=(samples, bands) del wvc, vis, vza, raa # Detect smile effects at atmospheric absorption features logger.info('Detect smile effects at atmosphere absorption features.') shifts = [] band_indices = [] n_features = 0 for wave, wave_range in absorption_features.items(): # Get sensor band range sensor_wave0, sensor_band0 = get_closest_wave(sensor_wave, wave_range[0]) sensor_wave1, sensor_band1 = get_closest_wave(sensor_wave, wave_range[1]) center_wave, band_index = get_closest_wave(sensor_wave, wave) # Check if continue if abs(sensor_wave0 - wave_range[0]) > 20 or abs(sensor_wave1 - wave_range[1]) > 20: continue # Get LUT band range _, atm_lut_band0 = get_closest_wave(atm_lut_wave, sensor_wave0 - 20) _, atm_lut_band1 = get_closest_wave(atm_lut_wave, sensor_wave1 + 20) # Optimize logger.info( 'Absorption feature center wavelength and range [nm] = %d: %d-%d.' % (wave, wave_range[0], wave_range[1])) x = [] for sample in range(samples): p = optimize.minimize( cost_fun, [0, 0], method='BFGS', args=(sensor_wave[sensor_band0:sensor_band1 + 1], sensor_fwhm[sensor_band0:sensor_band1 + 1], sensor_rdn[sensor_band0:sensor_band1 + 1, sample], atm_lut_wave[atm_lut_band0:atm_lut_band1], atm_lut_rdn[sample, atm_lut_band0:atm_lut_band1])) x.append(p.x) x = np.array(x) # Do linear interpolation for invalid values tmp_index = ~(np.abs(x[:, 0]) > 10.0) x[:, 0] = np.interp(np.arange(samples), np.arange(samples)[tmp_index], x[tmp_index, 0]) x[:, 1] = np.interp(np.arange(samples), np.arange(samples)[tmp_index], x[tmp_index, 1]) del tmp_index # Append shifts shifts.append(x) band_indices.append(band_index) n_features += 1 # Clear data del p, x, atm_lut_rdn sensor_rdn.flush() del sensor_rdn # Reshape data shifts = np.dstack(shifts).astype('float32').swapaxes(0, 1).swapaxes(1, 2) # Write shifts to binary file fid = open(sensor_dict['smile_effect_at_atm_features_file'], 'wb') fid.write(shifts[0, :, :].tostring()) # Shift in wavelength fid.write(shifts[1, :, :].tostring()) # Shift in fwhm fid.close() # Write header header = empty_envi_header() header[ 'description'] = 'Smile effects detected at atmosphere absorption features' header['file type'] = 'ENVI Standard' header['bands'] = 2 header['lines'] = n_features header['samples'] = samples header['interleave'] = 'bsq' header['byte order'] = 0 header['data type'] = 4 header['spectral center wavelengths'] = list(sensor_wave[band_indices]) header['spectral bandwiths'] = list(sensor_fwhm[band_indices]) header['band indices'] = band_indices write_envi_header( sensor_dict['smile_effect_at_atm_features_file'] + '.hdr', header) del header logger.info('Write smile effect at atmosphere absorption features to %s.' % sensor_dict['smile_effect_at_atm_features_file']) # Do interpolation fid = open(sensor_dict['smile_effect_file'], 'wb') x = np.arange(samples) y = sensor_wave[band_indices] x_new = x y_new = sensor_wave # Center wavelength z = np.zeros((shifts.shape[1], shifts.shape[2]), dtype='float64') for feature in range(shifts.shape[1]): p = np.poly1d(np.polyfit(x, shifts[0, feature, :], 4)) z[feature, :] = p(x) f = interpolate.interp2d(x, y, z, kind='cubic') z_new = f(x_new, y_new) + np.expand_dims(sensor_wave, axis=1) fid.write(z_new.astype('float32').tostring()) # Bandwidth z = np.zeros((shifts.shape[1], shifts.shape[2]), dtype='float64') for feature in range(shifts.shape[1]): p = np.poly1d(np.polyfit(x, shifts[1, feature, :], 5)) z[feature, :] = p(x) f = interpolate.interp2d(x, y, z, kind='cubic') z_new = f(x_new, y_new) + np.expand_dims(sensor_fwhm, axis=1) fid.write(z_new.astype('float32').tostring()) fid.close() del f, x, y, z, z_new # Write header header = empty_envi_header() header[ 'description'] = 'Spectral Center Wavelengths and Spectral Bandwidths' header['file type'] = 'ENVI Standard' header['bands'] = 2 header['lines'] = len(y_new) header['samples'] = len(x_new) header['interleave'] = 'bsq' header['byte order'] = 0 header['data type'] = 4 write_envi_header(sensor_dict['smile_effect_file'] + '.hdr', header) del header, x_new, y_new logger.info('Write smile effect to %s.' % sensor_dict['smile_effect_file'])
def make_quicklook(quicklook_figure_file, rdn_image_file, glt_image_file): """ Make a RGB quicklook image. Arguments: quicklook_figure_file: str Geo-rectified quicklook figure filename. rdn_image_file: str Radiance image filename, in BIL format. glt_image_file: str GLT image filename. """ if os.path.exists(quicklook_figure_file): logger.info('Save the quicklook figure to %s.' % quicklook_figure_file) return from ENVI import read_envi_header from Spectra import get_closest_wave # Read radiance image data. rdn_header = read_envi_header(os.path.splitext(rdn_image_file)[0] + '.hdr') rdn_image = np.memmap( rdn_image_file, # dtype='float32', dtype='int16', # to int, by Ting mode='r', shape=(rdn_header['lines'], rdn_header['bands'], rdn_header['samples'])) # Get RGB bands. if rdn_header['default bands'] is not None: rgb_bands = rdn_header['default bands'] else: rgb_bands = [] wave, _ = get_closest_wave(rdn_header['wavelength'], 450) if abs(wave - 450) < 10: for rgb_wave in [680, 550, 450]: _, band = get_closest_wave(rdn_header['wavelength'], rgb_wave) rgb_bands.append(band) else: for rgb_wave in [1220, 1656, 2146]: _, band = get_closest_wave(rdn_header['wavelength'], rgb_wave) rgb_bands.append(band) del band, wave # Read GLT image. glt_header = read_envi_header(os.path.splitext(glt_image_file)[0] + '.hdr') glt_image = np.memmap(glt_image_file, dtype=np.int32, mode='r', shape=(2, glt_header['lines'], glt_header['samples'])) # Write RGB image. driver = gdal.GetDriverByName('GTiff') ds = driver.Create(quicklook_figure_file, glt_header['samples'], glt_header['lines'], 3, gdal.GDT_Byte) ds.SetGeoTransform( (float(glt_header['map info'][3]), float(glt_header['map info'][5]), 0, float(glt_header['map info'][4]), 0, -float(glt_header['map info'][6]))) ds.SetProjection(glt_header['coordinate system string']) image = np.zeros((glt_header['lines'], glt_header['samples']), dtype='uint8') I, J = np.where((glt_image[0, :, :] >= 0) & (glt_image[1, :, :] >= 0)) for band_index, rgb_band in enumerate(rgb_bands): image[:, :] = 0 tmp_image = linear_percent_stretch(rdn_image[:, rgb_band, :]) image[I, J] = tmp_image[glt_image[0, I, J], glt_image[1, I, J]] ds.GetRasterBand(band_index + 1).WriteArray(image) del tmp_image glt_image.flush() rdn_image.flush() del glt_image, rdn_image ds = None del I, J, glt_header, rdn_header, image logger.info('Save the quicklook figure to %s.' % quicklook_figure_file)
def build_wvc_model(wvc_model_file, atm_lut_file, rdn_header_file, vis=40): """ Build water vapor models. Parameters ---------- wvc_model_file: str Water vapor column model filename. atm_lut_file: str Atmosphere lookup table file. rdn_header_file: str Radiance header filename. vis: float Visibility in [km]. """ if os.path.exists(wvc_model_file): logger.info('Write WVC models to %s.' % wvc_model_file) return from AtmLUT import read_binary_metadata, get_interp_range from ENVI import read_envi_header from Spectra import get_closest_wave, resample_spectra import json # Get sensor wavelengths & FWHMs header = read_envi_header(rdn_header_file) sensor_waves = np.array(header['wavelength']) sensor_fwhms = np.array(header['fwhm']) # Read atmospheric lookup table (ALT) metadata metadata = read_binary_metadata(atm_lut_file + '.meta') metadata['shape'] = tuple([int(v) for v in metadata['shape']]) atm_lut_WVC = np.array([float(v) for v in metadata['WVC']]) atm_lut_VIS = np.array([float(v) for v in metadata['VIS']]) atm_lut_VZA = np.array([float(v) for v in metadata['VZA']]) atm_lut_RAA = np.array([float(v) for v in metadata['RAA']]) atm_lut_WAVE = np.array([float(v) for v in metadata['WAVE']]) # VZA index vza_index = np.where(atm_lut_VZA == min(atm_lut_VZA))[0][-1] # RAA index raa_index = np.where(atm_lut_RAA == min(atm_lut_RAA))[0][-1] # Load atmospheric lookup table atm_lut = np.memmap( atm_lut_file, dtype='float32', mode='r', shape=metadata['shape']) # shape=(RHO, WVC, VIS, VZA, RAA, WAVE) # Subtract path radiance atm_lut_rdn = atm_lut[1, :, vza_index, raa_index, :] - atm_lut[ 0, :, vza_index, raa_index, :] # shape=(WVC, VIS, WAVE) atm_lut.flush() del atm_lut # VIS interpolation range vis_dict = get_interp_range(atm_lut_VIS, vis) # Do interpolation interp_rdn = np.zeros((len(atm_lut_WVC), len(atm_lut_WAVE))) for vis_index, vis_delta in vis_dict.items(): interp_rdn += atm_lut_rdn[:, vis_index, :] * vis_delta del atm_lut_rdn, vis_index, vis_delta # Build WVC models model_waves = [(890, 940, 1000), (1070, 1130, 1200)] wvc_model = dict() for waves in model_waves: # Get model wavelengths wave_1, wave_2, wave_3 = waves left_wave, left_band = get_closest_wave(sensor_waves, wave_1) middle_wave, middle_band = get_closest_wave(sensor_waves, wave_2) right_wave, right_band = get_closest_wave(sensor_waves, wave_3) # Build the model if np.abs(left_wave - wave_1) < 20 and np.abs( middle_wave - wave_2) < 20 and np.abs(right_wave - wave_3) < 20: model = dict() # Bands model['band'] = [int(left_band), int(middle_band), int(right_band)] # Wavelengths model['wavelength'] = [left_wave, middle_wave, right_wave] # Weights left_weight = (right_wave - middle_wave) / (right_wave - left_wave) right_weight = (middle_wave - left_wave) / (right_wave - left_wave) model['weight'] = [left_weight, right_weight] # Ratios & WVCs resampled_rdn = resample_spectra(interp_rdn, atm_lut_WAVE, sensor_waves[model['band']], sensor_fwhms[model['band']]) ratio = resampled_rdn[:, 1] / (left_weight * resampled_rdn[:, 0] + right_weight * resampled_rdn[:, 2]) index = np.argsort(ratio) model['ratio'] = list(ratio[index]) model['wvc'] = list(atm_lut_WVC[index]) # Save model parameters wvc_model['WVC_Model_%d' % wave_2] = model # Clear data del left_weight, right_weight, resampled_rdn, ratio, index, model del interp_rdn # Save WVC models to a .JSON file with open(wvc_model_file, 'w') as fid: json.dump(wvc_model, fid, indent=4) logger.info('Write WVC models to %s.' % wvc_model_file)
def estimate_wvc(wvc_file, rdn_file, sensors, sun_zenith, distance, background_mask_file): """ Estimate water vapor column. Parameters ---------- wvc_file: str Water vapor column image filename. rdn_file: str Radiance image filename. sensors: dict Sensors. mask_file: str Mask image filename. sun_zenith: float Sun zenith angle. distance: float Earth-to-sun distance. """ if os.path.exists(wvc_file): logger.info('Save the WVC image to %s.' % (wvc_file)) return from ENVI import read_envi_header, empty_envi_header, write_envi_header from Spectra import get_closest_wave, resample_solar_flux import json # Read radiance image rdn_header = read_envi_header(rdn_file + '.hdr') rdn_image = np.memmap(rdn_file, dtype='float32', mode='r', shape=(rdn_header['bands'], rdn_header['lines'], rdn_header['samples'])) # Read boundary mask bg_mask_header = read_envi_header(background_mask_file + '.hdr') bg_mask_image = np.memmap(background_mask_file, mode='r', dtype='bool', shape=(bg_mask_header['lines'], bg_mask_header['samples'])) # Mask out dark pixels good_mask_image = np.full((rdn_header['lines'], rdn_header['samples']), True) solar_flux = resample_solar_flux(solar_flux_file, rdn_header['wavelength'], rdn_header['fwhm']) cos_sun_zenith = np.cos(np.deg2rad(sun_zenith)) d2 = distance**2 # If reflectance at 850 nm is less than 0.01, mask out these pixels wave, band = get_closest_wave(rdn_header['wavelength'], 470) if abs(wave - 470) < 20: refl = rdn_image[band, :, :] * np.pi * d2 / (solar_flux[band] * cos_sun_zenith) good_mask_image &= (refl > 0.01) del refl # If reflectance at 850 nm is less than 0.10, mask out these pixels wave, band = get_closest_wave(rdn_header['wavelength'], 850) if abs(wave - 850) < 20: refl = rdn_image[band, :, :] * np.pi * d2 / (solar_flux[band] * cos_sun_zenith) good_mask_image &= (refl > 0.10) del refl # If reflectance at 1600 nm is less than 0.10, mask out these pixels wave, band = get_closest_wave(rdn_header['wavelength'], 1600) if abs(wave - 1600) < 20: refl = rdn_image[band, :, :] * np.pi * d2 / (solar_flux[band] * cos_sun_zenith) good_mask_image &= (refl > 0.10) del refl del wave, band, solar_flux, cos_sun_zenith, d2 # Estimate water vapor columns wvc_image = np.zeros((rdn_header['lines'], rdn_header['samples'])) for sensor_index, sensor_dict in sensors.items(): # Read water vapor column models wvc_model_file = sensor_dict['wvc_model_file'] wvc_models = json.load(open(wvc_model_file, 'r')) # Estimate WVC for wvc_model in wvc_models.values(): # Find band indices bands = [] for wave in wvc_model['wavelength']: _, band = get_closest_wave(rdn_header['wavelength'], wave) bands.append(band) del band # Calculate ratios ratio_image = rdn_image[bands[1], :, :] / ( 1e-10 + rdn_image[bands[0], :, :] * wvc_model['weight'][0] + rdn_image[bands[2], :, :] * wvc_model['weight'][1]) # Calculate water vapor columns wvc_image[good_mask_image] += np.interp( ratio_image[good_mask_image], wvc_model['ratio'], wvc_model['wvc']).astype('float32') # Clear data del ratio_image, bands, wvc_model rdn_image.flush() del rdn_image # Average wvc_image /= len(sensors) # Read solar flux wvc_image[~good_mask_image] = wvc_image[good_mask_image].mean() logger.info( 'WVC [mm] statistics: min=%.2f, max=%.2f, avg=%.2f, sd=%.2f.' % (wvc_image[good_mask_image].min(), wvc_image[good_mask_image].max(), wvc_image[good_mask_image].mean(), wvc_image[good_mask_image].std())) wvc_image[bg_mask_image] = -1000.0 del bg_mask_image, good_mask_image # Save water vapor column fid = open(wvc_file, 'wb') fid.write(wvc_image.astype('float32').tostring()) fid.close() # Write header wvc_header = empty_envi_header() wvc_header['description'] = 'Water vapor column [mm]' wvc_header['file type'] = 'ENVI Standard' wvc_header['bands'] = 1 wvc_header['lines'] = rdn_header['lines'] wvc_header['samples'] = rdn_header['samples'] wvc_header['file type'] = 'ENVI Standard' wvc_header['interleave'] = 'bsq' wvc_header['byte order'] = 0 wvc_header['data type'] = 4 wvc_header['data ignore value'] = -1000.0 wvc_header['map info'] = rdn_header['map info'] wvc_header['coordinate system string'] = rdn_header[ 'coordinate system string'] write_envi_header(wvc_file + '.hdr', wvc_header) del wvc_header, rdn_header logger.info('Save the WVC image to %s.' % (wvc_file))
def pre_classification(pre_class_image_file, rdn_image_file, sun_zenith, distance, background_mask_file=None): """ Pre-classify the image. Notes ----- (1) The classification algorithm used here is from ATCOR. Parameters ---------- pre_class_image_file: str Pre-classification image filename. rdn_image_file: str Radiance image filename, either BIL or BSQ. sun_zenith: float Sun zenith angle in degrees. distance: float Sun-Earth distance. background_mask_file: str Background mask filename. """ if os.path.exists(pre_class_image_file): logger.info('Write the pre-classification map to %s.' % (pre_class_image_file)) return from ENVI import empty_envi_header, read_envi_header, write_envi_header from Spectra import get_closest_wave, resample_solar_flux # Read radiance image data rdn_header = read_envi_header(os.path.splitext(rdn_image_file)[0] + '.hdr') if rdn_header['interleave'].lower() == 'bil': rdn_image = np.memmap(rdn_image_file, dtype='float32', mode='r', shape=(rdn_header['lines'], rdn_header['bands'], rdn_header['samples'])) elif rdn_header['interleave'].lower() == 'bsq': rdn_image = np.memmap(rdn_image_file, dtype='float32', mode='r', shape=(rdn_header['bands'], rdn_header['lines'], rdn_header['samples'])) else: logger.error('Cannot read radiance data from %s format file.' % rdn_header['interleave']) # Read solar flux solar_flux = resample_solar_flux(solar_flux_file, rdn_header['wavelength'], rdn_header['fwhm']) cos_sun_zenith = np.cos(np.deg2rad(sun_zenith)) d2 = distance**2 # Initialize classification pre_class_image = np.memmap(pre_class_image_file, dtype='uint8', mode='w+', shape=(rdn_header['lines'], rdn_header['samples'])) pre_class_image[:, :] = 0 # Define VNIR sensor wavelengths blue_wave, blue_band = get_closest_wave(rdn_header['wavelength'], 470) green_wave, green_band = get_closest_wave(rdn_header['wavelength'], 550) red_wave, red_band = get_closest_wave(rdn_header['wavelength'], 660) nir_wave, nir_band = get_closest_wave(rdn_header['wavelength'], 850) # Define SWIR sensor wavelengths cirrus_wave, cirrus_band = get_closest_wave(rdn_header['wavelength'], 1380) swir1_wave, swir1_band = get_closest_wave(rdn_header['wavelength'], 1650) swir2_wave, swir2_band = get_closest_wave(rdn_header['wavelength'], 2200) # Determine the sensor type if_vnir = abs(blue_wave - 470) < 20 and abs(green_wave - 550) < 20 and abs( red_wave - 660) < 20 and abs(nir_wave - 850) < 20 if_swir = abs(cirrus_wave - 1380) < 20 and abs( swir1_wave - 1650) < 20 and abs(swir2_wave - 2200) < 20 if_whole = if_vnir and if_swir # Do classification if if_whole: # Calculate the reflectance at different bands if rdn_header['interleave'].lower() == 'bil': blue_refl = rdn_image[:, blue_band, :] * np.pi * d2 / ( solar_flux[blue_band] * cos_sun_zenith) green_refl = rdn_image[:, green_band, :] * np.pi * d2 / ( solar_flux[green_band] * cos_sun_zenith) red_refl = rdn_image[:, red_band, :] * np.pi * d2 / ( solar_flux[red_band] * cos_sun_zenith) nir_refl = rdn_image[:, nir_band, :] * np.pi * d2 / ( solar_flux[nir_band] * cos_sun_zenith) cirrus_refl = rdn_image[:, cirrus_band, :] * np.pi * d2 / ( solar_flux[cirrus_band] * cos_sun_zenith) swir1_refl = rdn_image[:, swir1_band, :] * np.pi * d2 / ( solar_flux[swir1_band] * cos_sun_zenith) swir2_refl = rdn_image[:, swir2_band, :] * np.pi * d2 / ( solar_flux[swir2_band] * cos_sun_zenith) else: blue_refl = rdn_image[blue_band, :, :] * np.pi * d2 / ( solar_flux[blue_band] * cos_sun_zenith) green_refl = rdn_image[green_band, :, :] * np.pi * d2 / ( solar_flux[green_band] * cos_sun_zenith) red_refl = rdn_image[red_band, :, :] * np.pi * d2 / ( solar_flux[red_band] * cos_sun_zenith) nir_refl = rdn_image[nir_band, :, :] * np.pi * d2 / ( solar_flux[nir_band] * cos_sun_zenith) cirrus_refl = rdn_image[cirrus_band, :, :] * np.pi * d2 / ( solar_flux[cirrus_band] * cos_sun_zenith) swir1_refl = rdn_image[swir1_band, :, :] * np.pi * d2 / ( solar_flux[swir1_band] * cos_sun_zenith) swir2_refl = rdn_image[swir2_band, :, :] * np.pi * d2 / ( solar_flux[swir2_band] * cos_sun_zenith) rdn_image.flush() del rdn_image # Calculate NDSI ndsi = (green_refl - swir1_refl) / (green_refl + swir1_refl + 1e-10) # Water water = (nir_refl < 0.05) & (swir1_refl < 0.03) # Land land = ~water # Shadow shadow = (water | land) & (green_refl < 0.01) land = land & (~shadow) water = water & (~shadow) # Cloud over land Tc = 0.20 cloud_over_land = land & (blue_refl > Tc) & (red_refl > 0.15) & ( nir_refl < red_refl * 2.0) & (nir_refl > red_refl * 0.8) & ( nir_refl > swir1_refl) & (ndsi < 0.7) land = land & (~cloud_over_land) # Cloud over water cloud_over_water = water & (blue_refl < 0.40) & (blue_refl > 0.20) & ( green_refl < blue_refl) & (nir_refl < green_refl) & ( swir1_refl < 0.15) & (ndsi < 0.20) water = water & (~cloud_over_water) # Cloud shadow cloud_shadow = (water | land) & (nir_refl < 0.12) & ( nir_refl > 0.04) & (swir1_refl < 0.20) land = land & (~cloud_shadow) water = water & (~cloud_shadow) # Snow & ice snow_ice = land & (((blue_refl > 0.22) & (ndsi > 0.60)) | ((green_refl > 0.22) & (ndsi > 0.25) & (swir2_refl < 0.5 * green_refl))) land = land & (~snow_ice) # Cirrus over land & water thin_cirrus = (cirrus_refl < 0.015) & (cirrus_refl > 0.010) medium_cirrus = (cirrus_refl < 0.025) & (cirrus_refl >= 0.015) thick_cirrus = (cirrus_refl < 0.040) & (cirrus_refl >= 0.025) thin_cirrus_over_land = land & thin_cirrus land = land & (~thin_cirrus_over_land) medium_cirrus_over_land = land & medium_cirrus land = land & (~medium_cirrus_over_land) thick_cirrus_over_land = land & thick_cirrus land = land & (~thick_cirrus_over_land) thin_cirrus_over_water = water & thin_cirrus water = water & (~thin_cirrus_over_water) medium_cirrus_over_water = water & thin_cirrus water = water & (~medium_cirrus_over_water) thick_cirrus_over_water = water & thin_cirrus water = water & (~thick_cirrus_over_water) cirrus_cloud = (water | land) & (cirrus_refl < 0.050) & (cirrus_refl > 0.040) land = land & (~cirrus_cloud) water = water & (~cirrus_cloud) thick_cirrus_cloud = (water | land) & (cirrus_refl > 0.050) land = land & (~thick_cirrus_cloud) water = water & (~thick_cirrus_cloud) # Haze over water T2 = 0.04 T1 = 0.12 thin_haze_over_water = water & (nir_refl >= T1) & (nir_refl <= 0.06) water = water & (~thin_haze_over_water) medium_haze_over_water = water & (nir_refl >= 0.06) & (nir_refl <= T2) water = water & (~medium_haze_over_water) # Assign class values pre_class_image[shadow] = 1 pre_class_image[thin_cirrus_over_water] = 2 pre_class_image[medium_cirrus_over_water] = 3 pre_class_image[thick_cirrus_over_water] = 4 pre_class_image[land] = 5 pre_class_image[snow_ice] = 7 pre_class_image[thin_cirrus_over_land] = 8 pre_class_image[medium_cirrus_over_land] = 9 pre_class_image[thick_cirrus_over_land] = 10 pre_class_image[thin_haze_over_water] = 13 pre_class_image[medium_haze_over_water] = 14 pre_class_image[cloud_over_land] = 15 pre_class_image[cloud_over_water] = 16 pre_class_image[water] = 17 pre_class_image[cirrus_cloud] = 18 pre_class_image[thick_cirrus_cloud] = 19 pre_class_image[cloud_shadow] = 22 # Clear data del shadow, thin_cirrus_over_water, medium_cirrus_over_water del thick_cirrus_over_water, land, snow_ice del thin_cirrus_over_land, medium_cirrus_over_land, thick_cirrus_over_land del thin_haze_over_water, medium_haze_over_water, cloud_over_land del cloud_over_water, water, cirrus_cloud del thick_cirrus_cloud, cloud_shadow elif if_vnir: # Calculate the reflectance at different bands if rdn_header['interleave'].lower() == 'bil': blue_refl = rdn_image[:, blue_band, :] * np.pi * d2 / ( solar_flux[blue_band] * cos_sun_zenith) green_refl = rdn_image[:, green_band, :] * np.pi * d2 / ( solar_flux[green_band] * cos_sun_zenith) red_refl = rdn_image[:, red_band, :] * np.pi * d2 / ( solar_flux[red_band] * cos_sun_zenith) nir_refl = rdn_image[:, nir_band, :] * np.pi * d2 / ( solar_flux[nir_band] * cos_sun_zenith) else: blue_refl = rdn_image[blue_band, :, :] * np.pi * d2 / ( solar_flux[blue_band] * cos_sun_zenith) green_refl = rdn_image[green_band, :, :] * np.pi * d2 / ( solar_flux[green_band] * cos_sun_zenith) red_refl = rdn_image[red_band, :, :] * np.pi * d2 / ( solar_flux[red_band] * cos_sun_zenith) nir_refl = rdn_image[nir_band, :, :] * np.pi * d2 / ( solar_flux[nir_band] * cos_sun_zenith) rdn_image.flush() del rdn_image # Water water = nir_refl < 0.05 # Land land = ~water # Shadow shadow = (water | land) & (green_refl < 0.01) land = land & (~shadow) water = water & (~shadow) # Cloud over land Tc = 0.20 cloud_over_land = land & (blue_refl > Tc) & (red_refl > 0.15) & ( nir_refl < red_refl * 2.0) & (nir_refl > red_refl * 0.8) land = land & (~cloud_over_land) # Cloud over water cloud_over_water = water & (blue_refl < 0.40) & (blue_refl > 0.20) & ( green_refl < blue_refl) & (nir_refl < green_refl) water = water & (~cloud_over_water) # Cloud shadow cloud_shadow = (water | land) & (nir_refl < 0.12) & (nir_refl > 0.04) land = land & (~cloud_shadow) water = water & (~cloud_shadow) # Haze over water T2 = 0.04 T1 = 0.12 thin_haze_over_water = water & (nir_refl >= T1) & (nir_refl <= 0.06) water = water & (~thin_haze_over_water) medium_haze_over_water = water & (nir_refl >= 0.06) & (nir_refl <= T2) water = water & (~medium_haze_over_water) # Assign class values pre_class_image[shadow] = 1 pre_class_image[land] = 5 pre_class_image[thin_haze_over_water] = 13 pre_class_image[medium_haze_over_water] = 14 pre_class_image[cloud_over_land] = 15 pre_class_image[cloud_over_water] = 16 pre_class_image[water] = 17 pre_class_image[cloud_shadow] = 22 # Clear data del shadow, land, thin_haze_over_water del medium_haze_over_water, cloud_over_land, cloud_over_water del water, cloud_shadow elif if_swir: # Calculate the reflectance at different bands if rdn_header['interleave'] == 'bil': cirrus_refl = rdn_image[:, cirrus_band, :] * np.pi * d2 / ( solar_flux[cirrus_band] * cos_sun_zenith) swir1_refl = rdn_image[:, swir1_band, :] * np.pi * d2 / ( solar_flux[swir1_band] * cos_sun_zenith) swir2_refl = rdn_image[:, swir2_band, :] * np.pi * d2 / ( solar_flux[swir2_band] * cos_sun_zenith) else: cirrus_refl = rdn_image[cirrus_band, :, :] * np.pi * d2 / ( solar_flux[cirrus_band] * cos_sun_zenith) swir1_refl = rdn_image[swir1_band, :, :] * np.pi * d2 / ( solar_flux[swir1_band] * cos_sun_zenith) swir2_refl = rdn_image[swir2_band, :, :] * np.pi * d2 / ( solar_flux[swir2_band] * cos_sun_zenith) rdn_image.flush() del rdn_image # Water water = swir1_refl < 0.03 # Land land = ~water # Cirrus over land & water thin_cirrus = (cirrus_refl < 0.015) & (cirrus_refl > 0.010) medium_cirrus = (cirrus_refl < 0.025) & (cirrus_refl >= 0.015) thick_cirrus = (cirrus_refl < 0.040) & (cirrus_refl >= 0.025) thin_cirrus_over_land = land & thin_cirrus land = land & (~thin_cirrus_over_land) medium_cirrus_over_land = land & medium_cirrus land = land & (~medium_cirrus_over_land) thick_cirrus_over_land = land & thick_cirrus land = land & (~thick_cirrus_over_land) thin_cirrus_over_water = water & thin_cirrus water = water & (~thin_cirrus_over_water) medium_cirrus_over_water = water & thin_cirrus water = water & (~medium_cirrus_over_water) thick_cirrus_over_water = water & thin_cirrus water = water & (~thick_cirrus_over_water) cirrus_cloud = (water | land) & (cirrus_refl < 0.050) & (cirrus_refl > 0.040) land = land & (~cirrus_cloud) water = water & (~cirrus_cloud) thick_cirrus_cloud = (water | land) & (cirrus_refl > 0.050) land = land & (~thick_cirrus_cloud) water = water & (~thick_cirrus_cloud) # Assign class values pre_class_image[thin_cirrus_over_water] = 2 pre_class_image[medium_cirrus_over_water] = 3 pre_class_image[thick_cirrus_over_water] = 4 pre_class_image[land] = 5 pre_class_image[thin_cirrus_over_land] = 8 pre_class_image[medium_cirrus_over_land] = 9 pre_class_image[thick_cirrus_over_land] = 10 pre_class_image[water] = 17 pre_class_image[cirrus_cloud] = 18 pre_class_image[thick_cirrus_cloud] = 19 # Clear data del thin_cirrus_over_water, medium_cirrus_over_water, thick_cirrus_over_water, del land, thin_cirrus_over_land, medium_cirrus_over_land, thick_cirrus_over_land del water, cirrus_cloud, thick_cirrus_cloud else: logger.error( 'Cannot find appropriate wavelengths to do pre-classification.') # Apply background mask if background_mask_file is not None: bg_mask_header = read_envi_header( os.path.splitext(background_mask_file)[0] + '.hdr') bg_mask_image = np.memmap(background_mask_file, mode='r', dtype='bool', shape=(bg_mask_header['lines'], bg_mask_header['samples'])) pre_class_image[bg_mask_image] = 0 bg_mask_image.flush() del bg_mask_image # Clear data pre_class_image.flush() del pre_class_image # Write header pre_class_header = empty_envi_header() pre_class_header['description'] = 'Pre-classification' pre_class_header['file type'] = 'ENVI Standard' pre_class_header['samples'] = rdn_header['samples'] pre_class_header['lines'] = rdn_header['lines'] pre_class_header['classes'] = 23 pre_class_header['bands'] = 1 pre_class_header['byte order'] = 0 pre_class_header['header offset'] = 0 pre_class_header['interleave'] = 'bsq' pre_class_header['data type'] = 1 pre_class_header['class names'] = [ '0: geocoded background', '1: shadow', '2: thin cirrus (water)', '3: medium cirrus (water)', '4: thick cirrus (water)', '5: land (clear)', '6: saturated', '7: snow/ice', '8: thin cirrus (land)', '9: medium cirrus (land)', '10: thick cirrus (land)', '11: thin haze (land)', '12: medium haze (land)', '13: thin haze/glint (water)', '14: med. haze/glint (water)', '15: cloud (land)', '16: cloud (water)', '17: water', '18: cirrus cloud', '19: cirrus cloud (thick)', '20: bright', '21: topogr. shadow', '22: cloud shadow' ] pre_class_header['class lookup'] = [ 80, 80, 80, 0, 0, 0, 0, 160, 250, 0, 200, 250, 0, 240, 250, 180, 100, 40, 200, 0, 0, 250, 250, 250, 240, 240, 170, 225, 225, 150, 210, 210, 120, 250, 250, 100, 230, 230, 80, 80, 100, 250, 130, 130, 250, 180, 180, 180, 160, 160, 240, 0, 0, 250, 180, 180, 100, 160, 160, 100, 200, 200, 200, 40, 40, 40, 50, 50, 50 ] pre_class_header['map info'] = rdn_header['map info'] pre_class_header['coordinate system string'] = rdn_header[ 'coordinate system string'] write_envi_header( os.path.splitext(pre_class_image_file)[0] + '.hdr', pre_class_header) logger.info('Write the pre-classification map to %s.' % (pre_class_image_file))
def estimate_vis(vis_file, ddv_file, atm_lut_file, rdn_file, sca_file, background_mask_file): """ Estimate visibility. Arguments: vis_file: str Visibility map filename. ddv_file: str Dark dense vegetation map filename. atm_lut_file: str Atmospheric lookup table filename. rdn_file: str Radiance filename. sca_file: str Scan angle filename. background_mask_file: Background mask filename. """ if os.path.exists(ddv_file) and os.path.exists(vis_file): logger.info('Write the visibility map to %s.' % vis_file) logger.info('Write the DDV map to %s.' % ddv_file) return from ENVI import read_envi_header, empty_envi_header, write_envi_header from AtmLUT import read_binary_metadata from Spectra import get_closest_wave from AtmCorr import atm_corr_band # Read radiance header. rdn_header = read_envi_header(os.path.splitext(rdn_file)[0] + '.hdr') # Find VNIR and SWIR sensor wavelengths. red_wave, red_band = get_closest_wave(rdn_header['wavelength'], 660) nir_wave, nir_band = get_closest_wave(rdn_header['wavelength'], 850) swir1_wave, swir1_band = get_closest_wave(rdn_header['wavelength'], 1650) swir2_wave, swir2_band = get_closest_wave(rdn_header['wavelength'], 2130) # Determine the sensor type. if_vnir = abs(red_wave - 660) < 20 and abs(nir_wave - 850) < 20 if_swir = abs(swir1_wave - 1650) < 20 or abs(swir2_wave - 2130) < 20 if if_vnir and if_swir: logger.info( 'Both VNIR and SWIR bands are used for estimating visibility.') elif if_vnir: logger.info('Only VNIR bands are used for estimating visibility.') elif if_swir: logger.info( 'Only SWIR bands are available, a constant visibility (23 km) is used.' ) else: logger.error( 'Cannot find appropriate bands for estimating visibility.') # Read atmospheric lookup table. atm_lut_metadata = read_binary_metadata(atm_lut_file + '.meta') atm_lut_metadata['shape'] = tuple( [int(v) for v in atm_lut_metadata['shape']]) atm_lut_RHO = np.array([float(v) for v in atm_lut_metadata['RHO']]) atm_lut_WVC = np.array([float(v) for v in atm_lut_metadata['WVC']]) atm_lut_VIS = np.array([float(v) for v in atm_lut_metadata['VIS']]) atm_lut_VZA = np.array([float(v) for v in atm_lut_metadata['VZA']]) atm_lut_RAA = np.array([float(v) for v in atm_lut_metadata['RAA']]) atm_lut = np.memmap(atm_lut_file, dtype=atm_lut_metadata['dtype'], mode='r', shape=atm_lut_metadata['shape'] ) # shape = (RHO, WVC, VIS, VZA, RAA, WAVE) # Read radiance image. rdn_image = np.memmap( rdn_file, # dtype='float32', dtype='int16', # in int, by Ting mode='r', shape=(rdn_header['bands'], rdn_header['lines'], rdn_header['samples'])) # Read VZA and RAA image. sca_header = read_envi_header(os.path.splitext(sca_file)[0] + '.hdr') saa = float(sca_header['sun azimuth']) sca_image = np.memmap(sca_file, dtype='float32', offset=0, shape=(sca_header['bands'], sca_header['lines'], sca_header['samples'])) # vza vza_image = np.copy(sca_image[0, :, :]) # raa raa_image = saa - sca_image[1, :, :] raa_image[raa_image < 0] += 360.0 raa_image[raa_image > 180] = 360.0 - raa_image[raa_image > 180] raa_image[raa_image == 180] = 179.0 # Constrain RAA within bounds of ALT raa_max = atm_lut_RAA.max() raa_image[raa_image > raa_max] = raa_max - 5e-2 # clear data sca_image.flush() del sca_header, saa, sca_image # Set visibility and water vapor column values. metadata = read_binary_metadata(atm_lut_file + '.meta') tmp_wvc_image = np.ones( vza_image.shape) * atm_db_wvc_lut[metadata['atm_mode']] tmp_vis_image = np.ones(vza_image.shape) * 23 del metadata # Read background mask. bg_header = read_envi_header( os.path.splitext(background_mask_file)[0] + '.hdr') bg_mask = np.memmap(background_mask_file, dtype='bool', mode='r', shape=(bg_header['lines'], bg_header['samples'])) # Calculate NDVI. red_refl = atm_corr_band( atm_lut_WVC, atm_lut_VIS, atm_lut_VZA, atm_lut_RAA, atm_lut[..., red_band], tmp_wvc_image, tmp_vis_image, vza_image, raa_image, rdn_image[red_band, :, :] / 1000, # divided by 1000 to convert to original rdn, by Ting bg_mask) nir_refl = atm_corr_band(atm_lut_WVC, atm_lut_VIS, atm_lut_VZA, atm_lut_RAA, atm_lut[..., nir_band], tmp_wvc_image, tmp_vis_image, vza_image, raa_image, rdn_image[nir_band, :, :] / 1000, bg_mask) ndvi = (nir_refl - red_refl) / (nir_refl + red_refl + 1e-10) vis_image = np.zeros((rdn_header['lines'], rdn_header['samples'])) if if_vnir and if_swir: # Decide which SWIR band to use. if abs(swir2_wave - 2130) < 20: swir_wave = swir2_wave swir_band = swir2_band swir_refl_upper_bounds = [0.05, 0.10, 0.12] red_swir_ratio = 0.50 else: swir_wave = swir1_wave swir_band = swir1_band swir_refl_upper_bounds = [0.10, 0.15, 0.18] red_swir_ratio = 0.25 # Calculate swir refletance. swir_refl = atm_corr_band( atm_lut_WVC, atm_lut_VIS, atm_lut_VZA, atm_lut_RAA, atm_lut[..., swir_band], tmp_wvc_image, tmp_vis_image, vza_image, raa_image, rdn_image[swir_band, :, :] / 1000, # divided by 1000 to convert to original rdn, by Ting bg_mask) # Get DDV mask. for swir_refl_upper_bound in swir_refl_upper_bounds: ddv_mask = (ndvi > 0.10) & (swir_refl < swir_refl_upper_bound) & ( swir_refl > 0.01) percentage = np.sum(ddv_mask[~bg_mask]) / np.sum(~bg_mask) if percentage > 0.02: logger.info( 'The SWIR wavelength %.2f is used for detecting dark dense vegetation.' % swir_wave) logger.info('The SWIR reflectance upper boundary is %.2f.' % swir_refl_upper_bound) logger.info('The number of DDV pixels is %.2f%%.' % (percentage * 100)) break # Estimate Visibility. rows, columns = np.where(ddv_mask) # Estimate red reflectance. red_refl[rows, columns] = red_swir_ratio * swir_refl[rows, columns] del swir_refl for row, column in zip(rows, columns): interp_rdn = interp_atm_lut(atm_lut_RHO, atm_lut_WVC, atm_lut_VZA, atm_lut_RAA, atm_lut[..., red_band], red_refl[row, column], tmp_wvc_image[row, column], vza_image[row, column], raa_image[row, column]) vis_image[row, column] = np.interp( rdn_image[red_band, row, column] / 1000, interp_rdn[::-1], atm_lut_VIS[::-1] ) # divided by 1000 to convert to original rdn, by Ting rdn_image.flush() del rdn_image logger.info( 'Visibility [km] statistics: min=%.2f, max=%.2f, avg=%.2f, sd=%.2f.' % (vis_image[ddv_mask].min(), vis_image[ddv_mask].max(), vis_image[ddv_mask].mean(), vis_image[ddv_mask].std())) # Fill gaps with average values. vis_image[~ddv_mask] = vis_image[ddv_mask].mean() vis_image[bg_mask] = -1000.0 # Write the visibility data. fid = open(vis_file, 'wb') fid.write(vis_image.astype('float32').tostring()) fid.close() del vis_image vis_header = empty_envi_header() vis_header['description'] = 'Visibility [km]' vis_header['samples'] = rdn_header['samples'] vis_header['lines'] = rdn_header['lines'] vis_header['bands'] = 1 vis_header['byte order'] = 0 vis_header['header offset'] = 0 vis_header['interleave'] = 'bsq' vis_header['data ignore value'] = -1000.0 vis_header['data type'] = 4 vis_header['map info'] = rdn_header['map info'] vis_header['coordinate system string'] = rdn_header[ 'coordinate system string'] write_envi_header(os.path.splitext(vis_file)[0] + '.hdr', vis_header) logger.info('Write the visibility image to %s.' % vis_file) # Write the DDV data. fid = open(ddv_file, 'wb') fid.write(ddv_mask.tostring()) fid.close() del ddv_mask ddv_header = empty_envi_header() ddv_header[ 'description'] = 'DDV mask (1: Dark dense vegetaion; 0: non-dark dense vegetation)' ddv_header['samples'] = rdn_header['samples'] ddv_header['lines'] = rdn_header['lines'] ddv_header['bands'] = 1 ddv_header['byte order'] = 0 ddv_header['header offset'] = 0 ddv_header['interleave'] = 'bsq' ddv_header['data type'] = 1 ddv_header['map info'] = rdn_header['map info'] ddv_header['coordinate system string'] = rdn_header[ 'coordinate system string'] write_envi_header(os.path.splitext(ddv_file)[0] + '.hdr', ddv_header) logger.info('Write the DDV mask to %s.' % ddv_file) elif if_vnir: pass elif if_swir: pass # Clear data del ndvi, red_refl, nir_refl, rdn_header