Welcome to MIGHTI_L2’s documentation!

Contents:

MIGHTI_L2.attribute_measurement_location(tang_lat, tang_lon, tang_alt, integration_order=0)

Determine the geographical location to which the measurement will be attributed. Depending on integration_order (see function create_observation_matrix), this will either return the tangent locations, or the midpoint between two adjacent tangent locations.

INPUTS:

  • tang_lat – TYPE:array(ny), UNITS:deg. Tangent latitudes.
  • tang_lon – TYPE:array(ny), UNITS:deg. Tangent longitudes.
  • tang_alt – TYPE:array(ny), UNITS:km. Tangent altitudes.

OPTIONAL INPUTS:

  • integration_order – TYPE:int

    • 0: Use Riemann-sum rule for discretizing line-of-sight integral (default).
    • 1: Use trapezoidal rule for discretizing line-of-sight integral.

OUTPUTS:

  • lat – TYPE:array(ny), UNITS:deg. Measurement latitudes.
  • lon – TYPE:array(ny), UNITS:deg. Measurement longitudes.
  • alt – TYPE:array(ny), UNITS:km. Measurement altitudes.
MIGHTI_L2.bin_array(b, y, lon=False)

Downsample y by binning it, improving statistics. Every b elements of y will be averaged together to create a new array, y_b, of length ny_b = ceil(len(y)/b). Binning starts at the end of the array, so the first element of y_b may not represent exactly b samples of y.

INPUTS:

  • b – TYPE:int, The number of rows to bin together
  • y – TYPE:array(ny), The array to be binned

OPTIONAL INPUTS:

  • lon – TYPE:bool, If True, 360-deg discontinuities will

    be removed before averaging (e.g., for longitude binning).

OUTPUTS:

  • y_b – TYPE:array(ny_b), The binned array
MIGHTI_L2.bin_image(b, I)

Downsample the interferogram in altitude to improve statistics while degrading vertical resolution. Every b rows will be averaged together. Binning starts at high altitudes, so the lower rows of I_b may not represent exactly b rows of I.

INPUTS:

  • b – TYPE:int, The number of rows to bin together
  • I – TYPE:array(ny,nx), UNITS:arb. The MIGHTI interferogram

OUTPUTS:

  • I_b – TYPE:array(ny_b,nx), UNITS:arb. The binned MIGHTI interferogram
MIGHTI_L2.bin_uncertainty(b, ye)

Determine the uncertainty of a binned array from the uncertainty of the un-binned array. Specifically: If the array y has uncertainty given by array ye, then the array

y_b = bin_array(b, y)

has uncertainty given by the array

ye_b = bin_uncertainty(b, ye)

INPUTS:

  • b – TYPE:int, The number of rows to bin together
  • ye – TYPE:array(ny), The uncertainty of the pre-binned data

OUTPUTS:

  • ye_b – TYPE:array(ny_b), The uncertainty of the binned data
MIGHTI_L2.circular_mean(angle0, angle1)

Find the mean angle, taking into account 0/360 crossover. For example, circular_mean(10,50) is 30, but circular_mean(350,20) is 5.

INPUTS:

  • angle0 – TYPE:float or array, UNITS:deg. An angle in degrees.
  • angle1 – TYPE:float or array, UNITS:deg. An angle in degrees.

OUTPUTS:

  • angle – TYPE:float or array, UNITS:deg. The circular mean of the two input angles.
MIGHTI_L2.create_local_projection_matrix(tang_alt, icon_alt)

Define the matrix B whose entries give the factor by which a horizontal wind would be projected onto the line of sight. This has the same shape as the observation matrix (i.e., distance matrix). At the tangent point, this factor is 1.0. Far from the tangent point, this factor is smaller. If this effect is accounted for, it makes a small change in the winds (less than 5 m/s).

INPUTS:

  • tang_alt – TYPE:array(ny), UNITS:km. Tangent altitudes of each row of interferogram.
  • icon_alt – TYPE:float, UNITS:km. Altitude of the satellite.

OUTPUTS:

  • B – TYPE:array(ny,ny), UNITS:km. Local projection matrix. B[i,j] = cos(angle

    between ray i and the tangent of shell j at the point where they intersect)

MIGHTI_L2.create_observation_matrix(tang_alt, icon_alt, top_layer='exp', integration_order=0)

Define the matrix D whose inversion is known as “onion-peeling.”

The forward model is:

I = D * Ip

where I is the measured interferogram, D is the observation matrix, and Ip is the onion-peeled interferogram. If integration_order is 1, the observation matrix is created by assuming the spectrum (and thus the interferogram) is a piecewise linear function of altitude, treating the values of the interferogram at the tangent locations as the unknowns, and writing the measurements as a linear function of the unknowns. If integration_order is 0, the same recipe is followed, except the spectrum is assumed to be a piecewise constant function of altitude, and the unknowns are the values of the interferogram at the midpoint between two tangent altitudes.

Setting integration_order==0 is better for precision.

Setting integration_order==1 is better for accuracy.

INPUTS:

  • tang_alt – TYPE:array(ny), UNITS:km. Tangent altitudes of each row of interferogram.
  • icon_alt – TYPE:float, UNITS:km. Altitude of the satellite.

OPTIONAL INPUTS:

  • top_layer – TYPE:str, ‘thin’: assume VER goes to zero above top layer

    ‘exp’: assume VER falls off exponentially in altitude (default)

  • integration_order – TYPE:int,

    • 0: Use Riemann-sum rule for discretizing line-of-sight integral (default).
    • 1: Use trapezoidal rule for discretizing line-of-sight integral.

OUTPUTS:

  • D – TYPE:array(ny,ny), UNITS:km. Observation matrix. Also called the “path matrix”

    or “distance matrix”

MIGHTI_L2.extract_phase_from_row(row, zero_phase, phase_offset)

Given a 1-D interference pattern (i.e., a row of the intererogram), analyze it to get a single phase value, which represents the wind.

INPUTS:

  • row – TYPE:array(nx), UNITS:arb. A row of the complex-valued, MIGHTI interferogram.

  • zero_phase – TYPE:float, UNITS:rad. The phase angle which is equivalent

    to a wind value of zero.

  • phase_offset – TYPE:float, UNITS:rad. An offset to avoid 2pi ambiguities.

OUTPUTS:

  • phase – TYPE:float, UNITS:rad.
MIGHTI_L2.fix_longitudes(lons, lon_target)

Unwrap the list of longitudes to avoid 360-deg jumps. The list will be fixed so that it contains a value within 180 deg of lon_target and is otherwise continuous.

INPUTS:

  • lons – TYPE:array, UNITS:deg. An ordered list of longitudes to be unwrapped.
  • lon_target – TYPE:float, UNITS:deg. See above.

OUTPUTS:

  • lons_new – TYPE:array, UNITS:deg. An ordered list of longitudes with jumps removed.
MIGHTI_L2.interpolate_linear(x, y, x0, extrapolation='hold', prop_err=False, yerr=None)

Linear interpolation of the function y = f(x) to the location x0. x and y are vectors comprising samples of this function. There is also an option to propagate errors to the interpolated value. This function is 5 times faster than scipy.interpolate.interp1d, and allows for zero-order-hold extrapolation. If you are interpolating to many points, then scipy.interpolate.interp1d is likely faster.

INPUTS:

  • x – TYPE:array(n), UNITS:arb. Independent variable of samples of function.
  • y – TYPE:array(n), UNITS:arb. Dependent variable of samples of function.
  • x0 – TYPE:float, UNITS:arb. Independent variable of interpolation point.

OPTIONAL INPUTS:

  • extrapolation – TYPE:str, ‘hold’: extrapolate by using values at end points (default)

    ‘none’: do not extrapolate. Points will be np.nan

  • prop_err – TYPE:bool,

    • True: propagate errors from original to interpolated

      value, and return an extra output; yerr must be specified as an input.

    • False: do not propagate errors, and return only one

      output (default).

  • yerr – TYPE:array(n), UNITS:arb. Error in y, to be propagated to interpolated value.

OUTPUTS:

  • y0 – TYPE:float, UNITS:arb. Interpolated value.

OPTIONAL OUTPUT (if prop_err = True):

  • y0err – TYPE:float, UNTIS:arb. Propagated error of y0.
MIGHTI_L2.level1_dict_to_level21_dict(L1_dict, sigma, zero_phase, phase_offset, top_layer='exp', integration_order=0, account_for_local_projection=True, bin_size=1)

High-level function to run the Level 2.1 processing. It takes a dictionary (containing input variables extracted from a Level 1 file) and outputs a dictionary (containing output variables, which can be written to a file with another function).

INPUTS:

  • L1_dict – TYPE:dict. A dictionary containing variables needed for

    the level 2.1 processing:

    • L1_fn
    • I_amp
    • I_phase
    • I_amp_uncertainty
    • I_phase_uncertainty
    • tang_alt_start
    • tang_alt_stop
    • tang_lat_start
    • tang_lat_stop
    • tang_lon_start
    • tang_lon_stop
    • emission_color
    • icon_alt_start
    • icon_alt_stop
    • icon_lat_start
    • icon_lat_stop
    • icon_lon_start
    • icon_lon_stop
    • mighti_ecef_vectors_start
    • mighti_ecef_vectors_stop
    • icon_ecef_ram_vector_start
    • icon_ecef_ram_vector_stop
    • icon_velocity_start
    • icon_velocity_stop
    • source_files
    • time_start
    • time_stop
    • exp_time
    • optical_path_difference
  • sigma – TYPE:float, UNITS:m^-1. The wavenumber of the emission (1/wavelength)

  • zero_phase – TYPE:float, UNITS:rad. The phase corresponding to zero wind.

  • phase_offset – TYPE:float, UNITS.rad. A phase that is added to avoid 0/2pi discontinuities.

OPTIONAL INPUTS:

  • top_layer – TYPE:str, ‘thin’: assume VER goes to zero above top layer

    ‘exp’: assume VER falls off exponentially in altitude (default)

  • integration_order – TYPE:int, 0: Use Riemann-sum rule for discretizing line-of-sight integral (default)

    1: Use trapezoidal rule for discretizing line-of-sight integral

  • account_for_local_projection – TYPE:bool. If False, a simple inversion is used.

    If True, the inversion accounts for the fact that the ray is not perfectly tangent to each shell at each point along the ray. (default True)

  • bin_size – TYPE:int, The number of rows of the interferogram to bin together to

    improve statistics at the cost of altitude resolution. (default 1)

OUTPUTS:

  • L21_dict – TYPE:dict. A dictionary containing output variables of the Level 2.1 processing:

    • los_wind
    • los_wind_error
    • lat
    • lon
    • alt
    • time_start
    • time_stop
    • exp_time
    • az
    • emission_color
    • resolution_along_track
    • resolution_cross_track
    • resolution_alt
    • icon_alt
    • icon_lat
    • icon_lon
    • fringe_amplitude
    • fringe_amplitude_error
    • mighti_ecef_vectors
    • icon_velocity_ecef_vector
    • file_creation_time
    • source_files
    • bin_size
    • top_layer
    • integration_order
    • zero_phase
    • phase_offset
    • param_fn

TODO:

  • error handling
MIGHTI_L2.level1_to_dict(L1_fn, emission_color)

Read a level 1 file and translate it into a dictionary that the level 2.1 processing can use.

INPUTS:

  • L1_fn – TYPE:str. The full path and filename of the level 1 file.
  • emission_color – TYPE:str, ‘green’ or ‘red’.

OUTPUTS:

  • L1_dict – TYPE:dict. A dictionary containing information needed for

    the level 2.1 processing. See documentation for level1_dict_to_level21(...) for required keys.

TODO:

  • lots of stuff, mostly just waiting on finalization of L1 file.
  • time in better format
  • exposure time different than stop minus start?
  • will ICON velocity be specified in m/s eventually?
MIGHTI_L2.level1_to_level21(L1_fn, param_fn, emission_color, L21_path)

Highest-level function to apply the Level-1-to-Level-2.1 algorithm to a Level 1 file. The parameters of the algorithm are specified by a parameter file (in NetCDF4 format) with variables described below. There is a distinct parameter file for each sensor (A or B) and each emission color (red or green), but the algorithm is identical. The parameter file can be created with the “save_nc_param_file” function below.

INPUTS:

  • L1_fn – TYPE:str. The full path to the Level 1 file to be processed

  • param_fn – TYPE:str. The full path to the Level 2 parameters NetCDF4 file, with the following

    variables, described in more detail in the “save_nc_param_file” function below:

    • ICON_MIGHTI_Wavenumber
    • ICON_MIGHTI_Zero_Wind_Phase
    • ICON_MIGHTI_Phase_Offset
    • ICON_MIGHTI_Top_Layer
    • ICON_MIGHTI_Integration_Order
    • ICON_MIGHTI_Bin_Size
  • emission_color – TYPE:str. ‘red’ or ‘green’

  • L21_path – TYPE:str. The directory the Level 2.1 file will be saved in, including trailing “/”

    (e.g., ‘/home/user/’)

OUTPUTS:

  • L21_fn – TYPE:str. The full path to the saved L2.1 file.
MIGHTI_L2.level1_to_level21_without_param_file(L1_fn, emission_color, L21_path, sigma, zero_phase, phase_offset, top_layer='exp', integration_order=0, account_for_local_projection=True, bin_size=1, param_fn=None)

High-level function to apply the Level-1-to-Level-2.1 algorithm to a Level 1 file. This version of the function requires the user to input the parameters instead of specifying them via a parameter file. A single Level 1 file can become two Level 2.1 files (one per emission color).

INPUTS:

  • L1_fn – TYPE:str. The full path to the Level 1 file to be processed

  • emission_color – TYPE:str. ‘red’ or ‘green’

  • L21_path – TYPE:str. The directory the Level 2.1 file will be saved in, including trailing “/”

    (e.g., ‘/home/user/’)

  • sigma – TYPE:float, UNITS:m^-1. The wavenumber of the emission (1/wavelength)

  • zero_phase – TYPE:float, UNITS:rad. The phase corresponding to zero wind.

  • phase_offset – TYPE:float, UNITS.rad. A phase that is added to avoid 0/2pi discontinuities.

OPTIONAL INPUTS:

  • top_layer – TYPE:str, ‘thin’: assume VER goes to zero above top layer

    ‘exp’: assume VER falls off exponentially in altitude (default)

  • integration_order – TYPE:int, 0: Use Riemann-sum rule for discretizing line-of-sight integral (default)

    1: Use trapezoidal rule for discretizing line-of-sight integral

  • account_for_local_projection – TYPE:bool. If False, a simple inversion is used.

    If True, the inversion accounts for the fact that the ray is not perfectly tangent to each shell at each point along the ray. (default True)

  • bin_size – TYPE:int, The number of rows of the interferogram to bin together to

    improve statistics at the cost of altitude resolution. (default 1)

  • param_fn – TYPE:str, (default None). If specified, this string will be written to the

    L2.1 NetCDF4 file’s Calibration_File attribute. If this function is being called directly, the intention is to leave as None.

OUTPUTS:

  • L21_fn – TYPE:str. The full path to the saved L2.1 file.
MIGHTI_L2.level1_uiuc_to_dict(L1_uiuc_fn)

THIS MAY BE DEPRECATED

Read a level 1 file in special “UIUC” format, and translate it into a dictionary that the level 2.1 processing can use. This function is intended to be temporary, and used for testing only, not in operations. The L1_uiuc file is a .npz file.

INPUTS:

  • L1_uiuc_fn – TYPE:str. The full path and filename of the level 1 file

OUTPUTS:

  • L1_dict – TYPE:dict. A dictionary containing information needed for

    the level 2.1 processing. See documentation for level1_dict_to_level21(...) for required keys.

MIGHTI_L2.level21_dict_to_level22_dict(L21_A_dict, L21_B_dict, sph_asym_thresh)

Given Level 2.1 data from MIGHTI A and MIGHTI B, process it with the Level 2.2 algorithm. This entails interpolating line-of-sight wind measurements from MIGHTI A and B to a common grid, and rotating to a geographic coordinate system to derive vector horizontal winds.

INPUTS:

  • L21_A_dict – TYPE:dict. The dictionary corresponding to three orbits of MIGHTI A measurements

    for a single emission color. See “level21_to_dict” for required keys.

  • L21_B_dict – TYPE:dict. The dictionary corresponding to three orbits of MIGHTI B measurements

    for a single emission color, which is the same as for the A measurements. See “level21_to_dict” for required keys.

  • sph_asym_thresh – TYPE:float. Relative difference in emission rate measured by A and B, beyond which

    the spherical-asymmetry flag will be raised. Technically, it should be “fringe amplitude” instead of “emission rate” due to the temperature dependence. Relative difference is defined as abs(A-B)/((A+B)/2).

OUTPUTS:

  • L22_dict – TYPE:dict. A dictionary containing the following variables. Most are given as

    arrays with shape (ny,nx) where ny is the number of altitude grid points, and nx is the number of horizontal grid points.

    • lat – TYPE:array(ny,nx), UNITS:deg. Latitude of each point on the grid.

    • lon – TYPE:array(ny,nx), UNITS:deg. Longitude of each point on the grid.

    • alt – TYPE:array(ny,nx), UNITS:km. Altitude of each point on the grid.

    • u – TYPE:array(ny,nx), UNITS:m/s. Estimated zonal wind (positive eastward).

    • v – TYPE:array(ny,nx), UNITS:m/s. Estimated meridional wind (positive northward).

    • u_error – TYPE:array(ny,nx), UNITS:m/s. Uncertainty in u.

    • v_error – TYPE:array(ny,nx), UNITS:m/s. Uncertainty in v.

    • error_flags – TYPE:array(ny,nx,ne), UNITS:none. The error flags (either 0 or 1) for each point

      in the grid. Each point has a number of error flags, which are set to 1 under the following circumstances:

      • 0 = missing A file

      • 1 = missing B file

      • 2 = A signal too weak

      • 3 = B signal too weak

      • 4 = A did not sample this altitude

      • 5 = B did not sample this altitude

      • 6 = A sample exists but equals np.nan

      • 7 = B sample exists but equals np.nan

      • 8 = spherical asymmetry: A&B VER

        estimates disagree

      • 9 = unknown error

    • time – TYPE:array(ny,nx), UNITS:none. The average between the time of the MIGHTI A

      and B measurements that contribute to this grid point, given as a datetime object.

    • time_delta – TYPE:array(ny,nx), UNITS:s. The difference between the time of the MIGHTI

      A and B measurements that contribute to this grid point.

    • ver_rel_diff – TYPE:array(ny,nx), UNITS:none. The difference between the fringe amplitude

      of the MIGHTI A and B measurements that contribute to this grid point, divided by the mean. When this is high, it indicates that spherical asymmetry may be a problem.

    • emission_color – TYPE:str, ‘red’ or ‘green’.

    • source_files – TYPE:list of str, All the files used to create the data product,

      including the full paths.

    • param_fn – TYPE:str, If one was used, the full path to the parameter

      file will be inserted here by a higher-level function. If not, it should be left as None.

TODO:

  • only save the parent files that were actually used.
  • explicitly pass in which files are previous, current, and next orbit
  • Will an orbit be -180 to 180 deg or 0 to 360? Who decides this?
MIGHTI_L2.level21_to_dict(L21_fns)

Load a series of Level 2.1 files and return relevant variables in a dictionary.

INPUTS:

  • L21_fns – TYPE:list of str. The paths to the Level 2.1 files to be loaded.

OUTPUTS:

  • L21_dict – TYPE:dict. A dictionary containing the following variables. Most

    are provided as arrays of shape (ny,nt), where ny is the number of altitude samples and nt is the number of time samples.

    • lat – TYPE:array(ny,nt), UNITS:deg. Sample latitudes.

    • lon – TYPE:array(ny,nt), UNITS:deg. Sample longitudes.

    • alt – TYPE:array(ny,nt), UNITS:km. Sample altitudes.

    • los_wind – TYPE:array(ny,nt), UNITS:m/s. Line-of-sight wind component towards MIGHTI.

    • los_wind_error – TYPE:array(ny,nt), UNITS:m/s. Error in los_wind variable.

    • local_az – TYPE:array(ny,nt), UNITS:deg. Azimuth angle of vector pointing from

      MIGHTI towards the sample location, at the sample location (deg E of N).

    • amp – TYPE:array(ny,nt), UNITS:arb. Fringe amplitude at sample locations

    • time – TYPE:array(nt). Array of datetime objects, one per file.

    • exp_time – TYPE:array(nt), UNITS:sec. Exposure time of each sample.

    • emission_color – TYPE:str, ‘red’ or ‘green’.

    • source_files – TYPE:list of str, A copy of the input.

MIGHTI_L2.level21_to_level22(L21_path, L22_path, param_fn)

Highest-level function to apply the Level-2.1-to-Level-2.2 algorithm to a series of Level 2.1 files (in the L21_path directory) and create a Level 2.2 file (in the L22_path directory). The parameters of the algorithm are specified by a parameters file (in NetCDF4 format) with variables described below. The parameters file can be created with the “save_nc_param_file” function below.

INPUTS:

  • L21_path – TYPE:str. The directory containing the L2.1 files to be processed, including

    trailing “/” (e.g., ‘/home/user/’). This should contain a series of MIGHTI A and B files, over three orbits (previous, current, and next).

  • L22_path – TYPE:str. The directory the L2.2 file will be saved in, including trailing “/”

    (e.g., ‘/home/user/’)

  • param_fn – TYPE:str. The full path to the Level 2 parameters NetCDF4 file, with the following

    variables, described in more detail in the “save_nc_param_file” function above:

    • ICON_MIGHTI_Spherical_Asymmetry_Threshold

OUTPUTS:

  • L22_fn – TYPE:str. The full path to the saved L2.2 file.
MIGHTI_L2.level21_to_level22_without_param_file(L21_path, L22_path, sph_asym_thresh, param_fn=None)

High-level function to apply the Level-2.1-to-Level-2.2 algorithm to a series of Level 2.1 files (in the L21_path directory) and create a Level 2.2 file (in the L22_path directory). This version of the function requires the user to input parameters manually, instead of specifying a parameter file.

INPUTS:

  • L21_path – TYPE:str. The directory containing the L2.1 files to be processed, including

    trailing “/” (e.g., ‘/home/user/’). This should contain a series of MIGHTI A and B files, over three orbits (previous, current, and next).

  • L22_path – TYPE:str. The directory the L2.2 file will be saved in, including trailing “/”

    (e.g., ‘/home/user/’)

  • sph_asym_thresh – TYPE:float. Relative difference in emission rate measured by A and B, beyond which

    the spherical-asymmetry flag will be raised. Technically, it should be “fringe amplitude” instead of “emission rate” due to the temperature dependence. Relative difference is defined as abs(A-B)/((A+B)/2).

OPTIONAL INPUTS:

  • param_fn – TYPE:str, (default None). If specified, this string will be written to the

    L2.2 NetCDF4 file’s Calibration_File attribute. If this function is being called directly, the intention is to leave as None.

OUTPUTS:

  • L22_fn – TYPE:str. The full path to the saved L2.2 file.
MIGHTI_L2.los_az_angle(sat_latlonalt, lat, lon, alt)

Calculate the azimuth angle of the line of sight, evaluated at the measurement location (lat, lon, alt). Assumes WGS84 Earth.

INPUTS:

  • sat_latlonalt – TYPE:array(3), UNITS:(deg,deg,km). Satellite location in WGS84.
  • lat – TYPE:array(ny), UNITS:deg. Measurement latitudes.
  • lon – TYPE:array(ny), UNITS:deg. Measurement longitudes.
  • alt – TYPE:array(ny), UNITS:km. Measurement altitudes.

OUTPUTS:

  • az – TYPE:array(ny), UNITS:deg. Azimuth angle of line of sight

    from the satellite to the measurement location, evaluated at the measurement location. Degrees East of North.

MIGHTI_L2.perform_inversion(I, tang_alt, icon_alt, I_phase_uncertainty, I_amp_uncertainty, zero_phase, phase_offset, top_layer='exp', integration_order=0, account_for_local_projection=True)

Perform the onion-peeling inversion on the interferogram to return a new interferogram, whose rows refer to specific altitudes. In effect, this function undoes the integration along the line of sight.

INPUTS:

  • I – TYPE:array(ny,nx), UNITS:arb. The complex-valued, MIGHTI interferogram.

  • tang_alt – TYPE:array(ny), UNITS:km. Tangent altitudes of each row of interferogram.

  • icon_alt – TYPE:float, UNITS:km. Altitude of the satellite.

  • I_phase_uncertainty – TYPE:array(ny), UNITS:rad. Uncertainty in the unwrapped, mean phase of each row of I.

    This is provided in L1 file.

  • I_amp_uncertainty – TYPE:array(ny), UNITS:arb. Uncertainty in the summed amplitude of each row of I.

    This is provided in L1 file.

  • zero_phase – TYPE:float, UNITS:rad. The phase angle which is equivalent

    to a wind value of zero.

  • phase_offset– TYPE:float, UNITS:rad. An offset to avoid 2pi ambiguities.

OPTIONAL INPUTS:

  • top_layer – TYPE:str, ‘thin’: assume VER goes to zero above top layer

    ‘exp’: assume VER falls off exponentially in altitude (default)

  • integration_order – TYPE:int,

    • 0: Use Riemann-sum rule for discretizing line-of-sight integral (default).
    • 1: Use trapezoidal rule for discretizing line-of-sight integral.
  • account_for_local_projection – TYPE:bool. If False, a simple inversion is used.

    If True, the inversion accounts for the fact that the ray is not perfectly tangent to each shell at each point along the ray. (default True)

OUTPUTS:

  • Ip – TYPE:array(ny,nx), UNITS:arb. The complex-valued, onion-peeled interferogram.
  • phase – TYPE:array(ny), UNITS:rad. The unwrapped, mean phase of each row of Ip.
  • amp – TYPE:array(ny), UNITS:arb. The amplitude of each row of Ip.
  • phase_uncertainty – TYPE:array(ny), UNITS:rad. The uncertainty of phase
  • amp_uncertainty – TYPE:array(ny), UNITS:rad. The uncertainty of amp
MIGHTI_L2.phase_to_wind_factor(sigma_opd)

Return the value f that satisfies w = f*p, where w is a wind change and p is a phase change. dphi = 2*pi*OPD*sigma*v/c (Eq 1 of Englert et al 2007 Appl Opt., and Eq 2 of Harding et al. 2016 SSR)

INPUTS:

  • sigma_opd – TYPE:float, UNITS:none. sigma times opd: Optical path difference measured

    in wavelengths. If analyzing an entire row at once, the mean OPD of that row should be used.

OUTPUTS:

  • f – TYPE:float, UNITS:m/s/rad. Phase to wind factor, described above.
MIGHTI_L2.remove_Earth_rotation(v_inertial, az, lat, lon, alt)

Transform wind measurement from inertial coordinates to a reference frame rotating with the Earth. This can be thought of as “removing Earth rotation from the line-of-sight measurement.”

INPUTS:

  • v_inertial – TYPE:array(ny), UNITS:m/s. Line-of-sight velocity in inertial

    coordinates, positive towards MIGHTI.

  • az – TYPE:array(ny), UNITS:deg. Azimuth angle of line of sight

    from the satellite to the measurement location, evaluated at the measurement location. Degrees East of North. See los_az_angle() above.

  • lat – TYPE:array(ny), UNITS:deg. Measurement latitudes.

  • lon – TYPE:array(ny), UNITS:deg. Measurement longitudes.

  • alt – TYPE:array(ny), UNITS:km. Measurement altitudes.

OUTPUTS:

  • v – TYPE:array(ny), UNITS:m/s. Line-of-sight velocity in Earth-fixed

    coordinates, positive towards MIGHTI.

MIGHTI_L2.remove_satellite_velocity(I, sat_latlonalt, sat_velocity, sat_velocity_vector, mighti_vectors, sigma_opd)

Modify the interferogram to remove the effect of satellite velocity upon the phase.

INPUTS:

  • I – TYPE:array(ny,nx), UNITS:arb. The MIGHTI interferogram.

  • sat_latlonalt – TYPE:array(3), UNITS:(deg,deg,km). Satellite location in WGS84.

  • sat_velocity – TYPE:float, UNITS:m/s. ICON velocity.

  • sat_velocity_vector – TYPE:array(3), UNITS:none. The unit ECEF vector describing

    the direction of ICON’s velocity vector.

  • mighti_vectors – TYPE:array(ny,nx,3), UNITS:none. mighti_vectors[i,j,:] is a unit 3-vector in ECEF

    coordinates defining the look direction of pixel (i,j).

  • sigma_opd – TYPE:array(nx), UNITS:none. The optical path difference (measured in wavelengths)

    for each column of the interferogram.

OUTPUTS:

  • I – TYPE:array(ny,nx), UNITS:arb. The MIGHTI interferogram, corrected

    for the effects of satellite motion on the phase.

MIGHTI_L2.save_nc_level21(path, L21_dict)

Take the output of the Level 2.1 processing and save it as a NetCDF4 file in the official format. NetCDF4 file conventions taken from “Science Operations Center Data Product Conventions” Rev 0.4, authored by Tori Fae., including the new filenaming requirement that will presumably be in Rev 0.5.

INPUTS:

  • path – TYPE:str. The directory the file will be saved in, including trailing “/”

    (e.g., ‘/home/user/’)

  • L21_dict – TYPE:dict. A dictionary containing output variables of the Level 2.1 processing.

    See documentation for level1_dict_to_level21_dict(...) for details.

OUTPUTS:

  • L21_fn – TYPE:str. The full path to the saved file.

TO-DO:

  • Maybe: Fill in more notes for each variable
  • accept L2.1 filename as argument, from which I will extract what should go in Data_Version
  • Have a second set of eyes look at descriptions of each variable
  • Should dimensions be labeled the same as variables? Altitude/Vector/Epoch. Should Depend_0 point to vars or dims?
MIGHTI_L2.save_nc_level22(path, L22_dict)

Take the output of the Level 2.2 processing and save it as a NetCDF4 file in the official format. NetCDF4 file conventions taken from “Science Operations Center Data Product Conventions” Rev 0.4, authored by Tori Fae., including the new filenaming requirement that will presumably be in Rev 0.5.

INPUTS:

  • path – TYPE:str. The directory the file will be saved in, including trailing “/”

    (e.g., ‘/home/user/’)

  • L22_dict – TYPE:dict. A dictionary containing output variables of the Level 2.2 processing.

    See documentation for level21_dict_to_level22_dict(...) for details.

OUTPUTS:

  • L22_fn – TYPE:str. The full path to the saved file.

TO-DO:

  • Maybe: Fill in more notes for each variable
  • accept L2.2 filename as argument, from which I will extract Data_Version
  • Have a second set of eyes look at descriptions of each variable
  • Should dimensions be labeled the same as variables? Altitude/Vector/Epoch. Should Depend_0 point to vars or dims?
  • Time Resolution attribute? 30-60 seconds? Varied? N/A?
  • Add VER/fringe-amplitude estimate?
  • So far I haven’t used Depend_0 or Depend_1 because there isn’t an array for altitude or longitude. It’s not a regular grid, in general. I’m not sure how NetCDF4/ISTP wants to handle that.
MIGHTI_L2.save_nc_param_file(path, param_dict, sensor, emission_color, source_files=[], t=None)

Save a file that contains inversion parameters for the MIGHTI L2.1 and possibly L2.2 processing (e.g., zero-wind phase). The filename format is:

ICON_MIGHTI-<sensor>_<emission_color>_L2PARAM_YYYY-MM-DD_HH.MM.SS.nc,

where the sensor (A or B), emission color (RED or GREEN), and, optionally, the date, are specified by the user. For example:

ICON_MIGHTI-A_RED_L2PARAM_2015-01-01_00.15.15.nc

The parameters to be saved are contained in the cal_dict input dictionary.

INPUTS:

  • path – TYPE:str. The directory the file will be saved in, including trailing “/”

    (e.g., ‘/home/user/’)

  • param_dict – TYPE:dict. A dictionary containing inversion parameters, which are listed below

    and explained in more detail in the code below.

    • sigma – UNITS:m^-1
    • zero_phase – UNITS:rad.
    • phase_offset – UNITS:rad.
    • bin_size
    • integration_order
    • top_layer
    • spher_asym_thresh
  • sensor – TYPE:str. ‘A’ or ‘B’

  • emission_color – TYPE:str. ‘red’ or ‘green’

OPTIONAL INPUTS:

  • source_files – TYPE:list of str. A list of the filenames that were used to generated this parameter

    file. This may point to a series of L2.1 files generated during a zero wind maneuver, or possibly L1 files. (Default [])

  • t – TYPE:datetime. A naive datetime.datetime object for writing to the filename.

    If None (default), the current time will be used.

OUTPUTS:

  • cal_fn – TYPE:str. The full path to the saved file.
MIGHTI_L2.tang_alt_to_ze(tang_alt, sat_alt, RE)

Return the zenith angle(s) of the look direction(s) with given tangent altitude(s). Uses spherical Earth approximation.

INPUTS:

  • tang_alt – TYPE:array or float, UNITS:km. Tangent altitude(s) of the ray(s)
  • sat_alt – TYPE:float, UNITS:km. Satellite altitude (sat_alt > tang_alt)
  • RE – TYPE:float, UNITS:km. Effective radius of Earth

OUTPUT:

  • ze – TYPE:array or float, UNITS:deg. Zenith angle(s) of the ray(s)
MIGHTI_L2.unwrap(x)

Unwrap a monotonically increasing phase signal to remove -2*pi jumps. This is very similar to np.unwrap, but only unwraps negative jumps.

INPUTS:

  • x – TYPE:array, UNITS:rad. Signal that has -2*pi jumps to remove

OUTPUTS:

  • xnew – TYPE:array, UNITS:rad. Copy of x with -2*pi jumps removed
MIGHTI_L2.ze_to_tang_alt(ze, sat_alt, RE)

Return the tangent altitude(s) of the look direction(s) with given zenith angle(s). Uses spherical Earth approximation.

INPUTS:

  • ze – TYPE:array or float, UNITS:deg. Zenith angle(s) of the ray(s)
  • sat_alt – TYPE:float, UNITS:km. Satellite altitude
  • RE – TYPE:float, UNITS:km. Effective radius of Earth

OUTPUT:

  • tang_alt – TYPE:array or float, UNITS:km. Tangent altitude(s) of the ray(s)
ICON.azze_to_ecef(latlonalt, az, ze)

Convert a line of sight in (az,ze) coordinates, defined at the location given in WGS84 coordinates [latitude, longitude, altitude], to a unit vector in earth-centered earth-fixed (ECEF) coordinates [x,y,z]. This function is very similiar to ven_to_ecef.

INPUTS:
  • latlonalt - a length-3 array containing the WGS-84 coordinates of the location:

    [latitude (degrees), longitude (degrees), altitude (km)]

  • az - direction of line of sight. degrees east of north.

  • ze - direction of line of sight. degrees down from zenith

OUTPUTS:
  • xyz - the unit direction vector in ECEF coordinates [x,y,z]
HISTORY:
ICON.azze_to_lla(satlatlonalt, az, ze, ht)

Find the location (lat, lon, alt) where the requested look direction (az, ze) from the requested location intersects the specified altitude

INPUTS:
  • satlatlonalt - array [latitude (deg), longitude (deg), altitude (km)] of the satellite
  • az - direction of line of sight. degrees east of north.
  • ze - direction of line of sight. degrees down from zenith.
  • ht - altitude (km) to calculate the intersection point for.
OUTPUTS:
  • latlonalt - array [latitude (deg), longitude (deg), altitude (km)] of the tangent location
HISTORY:
  • 22-Jan-2015: Written by Jonathan J. Makela (jmakela@illinois.edu) based on ICON.tangent_point
  • 13-Oct-2015: Updated by Brian J. Harding to use distance_to_shell function.
ICON.distance_to_shell(satlatlonalt, az, ze, shell_altitude, intersection='first', tol=1e-05)

Along the line of sight, find the distance from the satellite to the first intersection with the shell at a given altitude. If no intersection exists, return np.nan. If the satellite is below this shell, return np.nan.

INPUTS:
  • satlatlonalt - array [latitude (deg), longitude (deg), altitude (km)] of the satellite
  • az - direction of line of sight. degrees east of north.
  • ze - direction of line of sight. degrees down from zenith
  • shell_altitude - the altitude of the shell in km.
OPTIONAL INPUTS:
  • intersection - str, ‘first’ or ‘second’. In general there are two intersections of the

    ray with the shell. This argument specifies which will be returned.

  • tol - km, stopping criterion for the solver. Stop when the solution is not changing

    by more than tol.

OUTPUTS:
  • d - distance in km. If no intersection exists, d is nan.
HISTORY:
  • 10-Dec-2014: Written by Brian J. Harding (bhardin2@illinois.edu)
  • 31-Mar-2015: Included “intersection” parameter.
ICON.distance_to_tangent_point(satlatlonalt, az, ze)

Along the line of sight, find the distance from the satellite to the tangent point.

INPUTS:
  • satlatlonalt - array [latitude (deg), longitude (deg), altitude (km)] of the satellite
  • az - direction of line of sight. degrees east of north.
  • ze - direction of line of sight. degrees down from zenith
OUTPUTS:
  • d - distance in km.
HISTORY:
ICON.earth_curvature(lat)

The radius of curvature of the Earth at the location specified. This is different from the radius of Earth at that point. WGS84 Earth is assumed. INPUTS:

  • lat – TYPE:float, UNITS:deg. Latitude of point on surface.

OUTPUTS:

  • r – TYPE:float, UNITS:km. Radius of curvature of Earth.
ICON.ecef_to_azze(latlonalt, ecef)

Convert a direction vector (e.g., velocity) in earth-centered earth-fixed (ECEF) coordinates [dx,dy,dz], defined at the location given in WGS84 coordinates [latitude, longitude, altitude], to the azimuth and zenith angles of the direction vector. This function is similar to ecef_to_ven.

INPUTS:
  • latlonalt - a length-3 array containing the WGS-84 coordinates of the location:

    [latitude (degrees), longitude (degrees), altitude (km)]

  • ecef - the direction vector [dx,dy,dz] which will be converted to

    azimuth and zenith angles.

OUTPUTS:
  • az,ze - the azimuth (degrees East from North) and zenith (degrees down

    from Vertical) angles of the direction vector.

HISTORY:
ICON.ecef_to_ven(latlonalt, ecef)

Convert a direction vector (e.g., velocity) in earth-centered earth-fixed (ECEF) coordinates [dx,dy,dz], defined at the location given in WGS84 coordinates [latitude, longitude, altitude], to a vector in local vertical-east-north (VEN) coordinates.

INPUTS:
  • latlonalt - a length-3 array containing the WGS-84 coordinates of the location:

    [latitude (degrees), longitude (degrees), altitude (km)]

  • ecef - the direction vector [dx,dy,dz] which will be converted to

    VEN coordinates.

OUTPUTS:
  • ven - the direction vector in [vertical, east, north] coordinates.
HISTORY:
  • 06-Jan-2015: Written by Brian J. Harding (bhardin2@illinois.edu)
  • 02-Sep-2015: Changed from numerical to analytical solution (BJH)
ICON.ecef_to_wgs84(ecef_xyz)

Convert from earth-centered earth-fixed (ECEF) coordinates (x,y,z) to WGS-84 latitude, longitude, and altitude. INPUTS:

  • ecef_xyz - a length-3 array containing the X, Y, and Z locations in ECEF

    coordinates in kilometers.

OUTPUTS:
  • latlonalt - a length-3 array containing the WGS-84 coordinates:

    [latitude (degrees), longitude (degrees), altitude (km)] Altitude is defined as the height above the reference ellipsoid. Longitude is defined in [0,360).

HISTORY:
ICON.get_solar_zenith_angle(pt)

Calculate the angle from zenith to the sun at the time and location of the pyglow.Point pt. INPUT:

  • pt - pyglow.Point

OUTPUT:

  • sza - Solar zenith angle (deg)

HISTORY:

  • 13-Oct-2015: Written by Brian J. Harding
ICON.project_line_of_sight(satlatlonalt, az, ze, step_size=1.0, total_distance=4000.0)

Starting at the satellite, step along a line of sight, and return an array of equally-spaced points at intervals along this line of sight, in WGS84 latitude, longitude, and altitude coordinates.

INPUTS:
  • satlatlonalt - array [latitude (deg), longitude (deg), altitude (km)] of the satellite
  • az - direction of line of sight. degrees east of north.
  • ze - direction of line of sight. degrees down from zenith
OPTIONAL INPUTS:
  • step_size - spacing between points in the returned array (km).
  • total_distance - length of the projected line (km).
OUTPUTS:
  • xyz - array of size 3xN, where N is floor(total_distance/step_size). Contains

    the ECEF coordinates of every point along the line of sight, in step_size intervals.

  • latlonalt - array of size 3xN, where each column contains the latitude, longitude

    and altitude corresponding to the column of xyz, in WGS84 coordinates

HISTORY:
ICON.tangent_point(satlatlonalt, az, ze, tol=1e-07)

Find the location (lat, lon, alt) of the tangent point of a ray from the satellite. Current implementation finds the minimum altitude along the line of sight by using a numerical method with a precision of about 0.1 m. The vertical precision is much better than the horizontal precision.

INPUTS:
  • satlatlonalt - array [latitude (deg), longitude (deg), altitude (km)] of the satellite
  • az - direction of line of sight. degrees east of north.
  • ze - direction of line of sight. degrees down from zenith
OPTIONAL INPUTS:
  • tol - km, stopping criterion for the solver. Stop when the solution is not moving

    by more than a horizontal distance of tol.

OUTPUTS:
  • latlonalt - array [latitude (deg), longitude (deg), altitude (km)] of the tangent location.
HISTORY:
  • 10-Dec-2014: Written by Brian J. Harding (bhardin2@illinois.edu)
  • 02-Sep-2015: Use scipy.optimize to do minimization instead of my own implementation.
ICON.ven_to_ecef(latlonalt, ven)

Convert a direction vector (e.g., velocity) in local vertical-east-north (VEN) coordinates, defined at the location given in WGS84 coordinates [latitude, longitude, altitude], to a vector in earth-centered earth-fixed (ECEF) coordinates [x,y,z].

INPUTS:
  • latlonalt - a length-3 array containing the WGS-84 coordinates of the location:

    [latitude (degrees), longitude (degrees), altitude (km)]

  • ven - the local direction vector [vertical, east, north] which will be converted to

    ECEF coordinates.

OUTPUTS:
  • xyz - the direction vector in ECEF coordinates [x,y,z]
HISTORY:
  • 10-Nov-2014: Written by Brian J. Harding (bhardin2@illinois.edu)
  • 02-Sep-2015: Changed from numerical to analytical solution (BJH)
ICON.wgs84_to_ecef(latlonalt)

Convert from WGS84 coordinates [latitude, longitude, altitude] to earth-centered earth-fixed coordinates (ECEF) [x,y,z]

INPUTS:
  • latlonalt - a length-3 array containing the WGS-84 coordinates:

    [latitude (degrees), longitude (degrees), altitude (km)]

OUTPUTS:
  • ecef_xyz - a length-3 array containing the X, Y, and Z locations in ECEF

    coordinates in kilometers.

HISTORY:
ICON.wgs84constants()

Return the constants associated with the WGS-84 coordinate system

OUTPUTS:

  • a - semi-major axis of Earth, km
  • b - semi-minor axis of Earth, km
  • e - eccentricity of Earth

Indices and tables