isofit.inversion.inverse_simple

heuristic_atmosphere(RT, instrument, x_RT, x_instrument, meas, geom)[source]

From a given radiance, estimate atmospheric state with band ratios. Used to initialize gradient descent inversions.

Parameters:
  • RT (RadiativeTransfer) – radiative transfer model to use

  • instrument (Instrument) – instrument for noise characterization

  • x_RT (array) – radiative transfer portion of the state vector

  • x_instrument (array) – instrument portion of the state vector

  • meas (array) – a one-D numpy vector of radiance in uW/nm/sr/cm2

  • geom (Geometry) – geometry object corresponding to given measurement

Returns:

updated estimate of x_RT

Return type:

x_new

invert_algebraic(surface, RT, instrument, x_surface, x_RT, x_instrument, meas, geom)[source]

Inverts radiance algebraically using Lambertian assumptions to get a reflectance.

Parameters:
  • surface (Surface) – surface model

  • RT (RadiativeTransfer) – radiative transfer model to use

  • instrument (Instrument) – instrument model

  • x_surface (array) – surface portion of the state vector

  • x_RT (array) – radiative transfer portion of the state vector

  • x_instrument (array) – instrument portion of the state vector

  • meas (array) – a one-D numpy vector of radiance in uW/nm/sr/cm2

  • geom (Geometry) – geometry object corresponding to given measurement

Returns:

estimate of the surface reflectance based on the given surface model and specified atmospheric state Ls: estimate of the emitted surface leaving radiance coeffs: atmospheric parameters for the forward model

Return type:

rfl_est

invert_analytical(fm, winidx, meas, geom, x_RT, num_iter=1, hash_table=None, hash_size=None, diag_uncert=True, outside_ret_const=-0.01)[source]

Perform an analytical estimate of the conditional MAP estimate for a fixed atmosphere. Based on the “Inner loop” from Susiluoto, 2022.

Parameters:
  • fm (ForwardModel) – isofit forward model

  • winidx (array) – indices of surface components of state vector (to be solved)

  • meas (array) – a one-D numpy vector of radiance in uW/nm/sr/cm2

  • geom (Geometry) – geometry object corresponding to given measurement

  • x_RT (array) – the radiative transfer state variables

  • num_iter (int) – number of interations to run through

  • hash_table (Optional[OrderedDict]) – a hash table to use locally

  • hash_size (Optional[int]) – max size of given hash table

  • diag_uncert (bool) – 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, meas, geom)[source]

Find an initial guess at the state vector. This currently uses traditional (non-iterative, heuristic) atmospheric correction.

Parameters:
  • forward (ForwardModel) – isofit forward model

  • meas (array) – a one-D numpy vector of radiance in uW/nm/sr/cm2

  • geom (Geometry) – 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, wl, l_shoulder=850, r_shoulder=1100, lw_init=(0.02, 0.3, 0.0002), lw_bounds=([0, 0.5], [0, 1.0], [-0.0004, 0.0004]), ewt_detection_limit=0.5, return_abs_co=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 (array) – surface reflectance spectrum

  • wl (array) – instrument wavelengths, must be same size as rfl_meas

  • l_shoulder (float) – wavelength of left absorption feature shoulder

  • r_shoulder (float) – wavelength of right absorption feature shoulder

  • lw_init (tuple) – initial guess for liquid water path length, intercept, and slope

  • lw_bounds (tuple) – lower and upper bounds for liquid water path length, intercept, and slope

  • ewt_detection_limit (float) – upper detection limit for ewt

  • return_abs_co (bool) – 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