isofit.inversion.inverse_simple =============================== .. py:module:: isofit.inversion.inverse_simple Functions --------- .. autoapisummary:: isofit.inversion.inverse_simple.heuristic_atmosphere isofit.inversion.inverse_simple.invert_algebraic isofit.inversion.inverse_simple.invert_analytical isofit.inversion.inverse_simple.invert_simple isofit.inversion.inverse_simple.invert_liquid_water isofit.inversion.inverse_simple.beer_lambert_model Module Contents --------------- .. py:function:: 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) From a given radiance, estimate atmospheric state with band ratios. Used to initialize gradient descent inversions. :param fm: isofit forward model :param x_surface: surface portion of the state vector :param x_RT: radiative transfer portion of the state vector :param x_instrument: instrument portion of the state vector :param meas: a one-D numpy vector of radiance in uW/nm/sr/cm2 :param geom: geometry object corresponding to given measurement :param wl_lo: Low wavelength to use for the continuum removal H2O fit :param wl_center: Center wavelength to use for the continuum removal H2O fit :param wl_hi: High wavelength to use for the continuum removal H2O fit :returns: updated estimate of x_RT :rtype: x_new .. py:function:: 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) 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. :param surface: surface model :param RT: radiative transfer model to use :param instrument: instrument model :param x_surface: surface portion of the state vector :param x_RT: radiative transfer portion of the state vector :param x_instrument: instrument portion of the state vector :param meas: a one-D numpy vector of radiance in uW/nm/sr/cm2 :param 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 :rtype: rfl_est .. py:function:: 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) Perform an analytical estimate of the conditional MAP estimate for a fixed atmosphere. Based on the "Inner loop" from Susiluoto, 2022. :param fm: isofit forward model :param winidx: indices of surface components of state vector (to be solved) :param meas: a one-D numpy vector of radiance in uW/nm/sr/cm2 :param geom: geometry object corresponding to given measurement :param x0: the initialization state including surface from the superpixel and the atm from the smoothed atmosphere. :param num_iter: number of interations to run through :param hash_table: a hash table to use locally :param hash_size: max size of given hash table :param diag_uncert: flag indicating whether to diagonalize the uncertainty :returns: MAP estimate of the mean S: diagonal conditional posterior covariance estimate :rtype: x .. py:function:: invert_simple(forward: ForwardModel, meas: numpy.array, geom: Geometry) Find an initial guess at the state vector. This currently uses traditional (non-iterative, heuristic) atmospheric correction. :param forward: isofit forward model :param meas: a one-D numpy vector of radiance in uW/nm/sr/cm2 :param geom: geometry object corresponding to given measurement :returns: estimate of the full statevector based on initial conditions, geometry, and a heuristic guess :rtype: x .. py:function:: 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) Given a reflectance estimate, fit a state vector including liquid water path length based on a simple Beer-Lambert surface model. :param rfl_meas: surface reflectance spectrum :param wl: instrument wavelengths, must be same size as rfl_meas and in units of nm :param l_shoulder: wavelength of left absorption feature shoulder :param r_shoulder: wavelength of right absorption feature shoulder :param lw_init: initial guess for liquid water path length, intercept, and slope :param lw_bounds: lower and upper bounds for liquid water path length, intercept, and slope :param ewt_detection_limit: upper detection limit for ewt :param 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 :rtype: solution .. py:function:: beer_lambert_model(x, y, wl, alpha_lw) 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. :param x: state vector (liquid water path length, intercept, slope) :param y: measurement (surface reflectance spectrum) :param wl: instrument wavelengths :param alpha_lw: wavelength dependent absorption coefficients of liquid water :returns: residual between modeled and measured surface reflectance :rtype: resid