isofit.inversion.inverse_simple
Functions
|
From a given radiance, estimate atmospheric state with band ratios. |
|
Inverts radiance algebraically using Lambertian assumptions to get a |
|
Perform an analytical estimate of the conditional MAP estimate for |
|
Find an initial guess at the state vector. This currently uses |
|
Given a reflectance estimate, fit a state vector including liquid water path length |
|
Function, which computes the vector of residuals between measured and modeled surface reflectance optimizing |
Module Contents
- heuristic_atmosphere(fm: ForwardModel, x_surface: numpy.array, x_RT: numpy.array, x_instrument: numpy.array, meas: numpy.array, geom: Geometry, wl_lo: int = 865, wl_center: int = 945, wl_hi: int = 1040)[source]
From a given radiance, estimate atmospheric state with band ratios. Used to initialize gradient descent inversions.
- Parameters:
fm – isofit forward model
x_surface – surface portion of the state vector
x_RT – radiative transfer portion of the state vector
x_instrument – instrument portion of the state vector
meas – a one-D numpy vector of radiance in uW/nm/sr/cm2
geom – geometry object corresponding to given measurement
wl_lo – Low wavelength to use for the continuum removal H2O fit
wl_center – Center wavelength to use for the continuum removal H2O fit
wl_hi – High wavelength to use for the continuum removal H2O fit
- Returns:
updated estimate of x_RT
- Return type:
x_new
- invert_algebraic(surface: Surface, RT: RadiativeTransfer, instrument: Instrument, x_surface: numpy.array, x_RT: numpy.array, x_instrument: numpy.array, meas: numpy.array, geom: Geometry)[source]
Inverts radiance algebraically using Lambertian assumptions to get a reflectance. If the appropriate transmittance terms are available, surface slope will enter for the initial guess. Otherwise, the surface will be treated as if flat.
- Parameters:
surface – surface model
RT – radiative transfer model to use
instrument – instrument model
x_surface – surface portion of the state vector
x_RT – radiative transfer portion of the state vector
x_instrument – instrument portion of the state vector
meas – a one-D numpy vector of radiance in uW/nm/sr/cm2
geom – geometry object corresponding to given measurement
- Returns:
estimate of the surface reflectance based on the given surface model and specified atmospheric state coeffs: atmospheric parameters used for the inversion, returned for convenience
- Return type:
rfl_est
- invert_analytical(fm: ForwardModel, winidx: numpy.array, meas: numpy.array, geom: Geometry, x0: numpy.array, sub_state, num_iter: int = 1, hash_table: OrderedDict = None, hash_size: int = None, diag_uncert: bool = True, outside_ret_const: float = -0.01, fill_value: float = -9999.0)[source]
Perform an analytical estimate of the conditional MAP estimate for a fixed atmosphere. Based on the “Inner loop” from Susiluoto, 2022.
- Parameters:
fm – isofit forward model
winidx – indices of surface components of state vector (to be solved)
meas – a one-D numpy vector of radiance in uW/nm/sr/cm2
geom – geometry object corresponding to given measurement
x0 – the initialization state including surface from the superpixel and the atm from the smoothed atmosphere.
num_iter – number of interations to run through
hash_table – a hash table to use locally
hash_size – max size of given hash table
diag_uncert – flag indicating whether to diagonalize the uncertainty
- Returns:
MAP estimate of the mean S: diagonal conditional posterior covariance estimate
- Return type:
x
- invert_simple(forward: ForwardModel, meas: numpy.array, geom: Geometry)[source]
Find an initial guess at the state vector. This currently uses traditional (non-iterative, heuristic) atmospheric correction.
- Parameters:
forward – isofit forward model
meas – a one-D numpy vector of radiance in uW/nm/sr/cm2
geom – geometry object corresponding to given measurement
- Returns:
estimate of the full statevector based on initial conditions, geometry, and a heuristic guess
- Return type:
x
- invert_liquid_water(rfl_meas: numpy.array, wl: numpy.array, l_shoulder: float = 850, r_shoulder: float = 1100, lw_init: tuple = (0.02, 0.3, 0.0002), lw_bounds: tuple = ([0, 0.5], [0, 1.0], [-0.0004, 0.0004]), ewt_detection_limit: float = 0.5, return_abs_co: bool = False)[source]
Given a reflectance estimate, fit a state vector including liquid water path length based on a simple Beer-Lambert surface model.
- Parameters:
rfl_meas – surface reflectance spectrum
wl – instrument wavelengths, must be same size as rfl_meas and in units of nm
l_shoulder – wavelength of left absorption feature shoulder
r_shoulder – wavelength of right absorption feature shoulder
lw_init – initial guess for liquid water path length, intercept, and slope
lw_bounds – lower and upper bounds for liquid water path length, intercept, and slope
ewt_detection_limit – upper detection limit for ewt
return_abs_co – if True, returns absorption coefficients of liquid water
- Returns:
estimated liquid water path length, intercept, and slope based on a given surface reflectance
- Return type:
solution
- beer_lambert_model(x, y, wl, alpha_lw)[source]
Function, which computes the vector of residuals between measured and modeled surface reflectance optimizing for path length of surface liquid water based on the Beer-Lambert attenuation law.
- Parameters:
x – state vector (liquid water path length, intercept, slope)
y – measurement (surface reflectance spectrum)
wl – instrument wavelengths
alpha_lw – wavelength dependent absorption coefficients of liquid water
- Returns:
residual between modeled and measured surface reflectance
- Return type:
resid