asp_plot.altimetry

Contents

asp_plot.altimetry#

Attributes#

Classes#

AlignmentResult

Outcome of Altimetry.align_and_evaluate().

Altimetry

Process and analyze ICESat-2 ATL06-SR altimetry data with ASP DEMs.

Functions#

gds_query_async(query_type, bounds, results_code[, email])

Submit an async query to the ODE GDS REST API.

Module Contents#

class asp_plot.altimetry.AlignmentResult#

Outcome of Altimetry.align_and_evaluate().

Plain dataclass — importing this or the method does not pull in any report/fpdf dependencies, so it is safe to use from notebooks.

status#
One of:
  • "insufficient_points": not enough ATL06-SR points for pc_align to run (the aligned DEM is removed if one was written).

  • "no_improvement": pc_align ran but p50 did not improve toward 0 by more than the improvement_threshold_pct; the aligned DEM has been removed.

  • "success": p50 improved by more than the threshold; the aligned DEM is retained and Altimetry.aligned_dem_fn points to it.

Type:

str

alignment_report_df#

The alignment report table produced by Altimetry.alignment_report(). Empty for insufficient_points.

Type:

pandas.DataFrame

aligned_dem_fn#

Path to the aligned DEM on success; None otherwise (the file is cleaned up on non-success branches).

Type:

str or None

improvement_pct#

(p50_beg - p50_end) / p50_beg * 100 when computable, else None.

Type:

float or None

message#

Short human-readable summary suitable for a status line in the PDF report.

Type:

str

parameters_used#

The kwargs used for the alignment attempt (e.g. processing_level, minimum_points, improvement_threshold_pct). Echoed into the report.

Type:

dict

aligned_dem_fn: str | None = None#
alignment_report_df: pandas.DataFrame | None = None#
improvement_pct: float | None = None#
message: str = ''#
parameters_used: dict#
status: str#
class asp_plot.altimetry.Altimetry(directory, dem_fn, aligned_dem_fn=None, atl06sr_processing_levels={}, atl06sr_processing_levels_filtered={}, **kwargs)#

Process and analyze ICESat-2 ATL06-SR altimetry data with ASP DEMs.

This class provides functionality to request, filter, and analyze ICESat-2 ATL06-SR altimetry data in conjunction with ASP-generated DEMs. It includes methods to request data from the SlideRule API, filter based on various criteria, align DEMs to altimetry data, and visualize the results.

directory#

Root directory for outputs and analysis

Type:

str

dem_fn#

Path to the DEM file to analyze

Type:

str

aligned_dem_fn#

Path to the aligned DEM file if available

Type:

str or None

atl06sr_processing_levels#

Dictionary of ATL06-SR data for different processing levels

Type:

dict

atl06sr_processing_levels_filtered#

Dictionary of filtered ATL06-SR data for different processing levels

Type:

dict

alignment_report_df#

DataFrame containing alignment reports when available

Type:

pandas.DataFrame or None

Examples

>>> altimetry = Altimetry('/path/to/directory', '/path/to/dem.tif')
>>> altimetry.request_atl06sr_multi_processing(save_to_parquet=True)
>>> altimetry.filter_esa_worldcover(filter_out="water")
>>> altimetry.alignment_report()
>>> altimetry.mapview_plot_atl06sr_to_dem()
align_and_evaluate(processing_level='all', improvement_threshold_pct=5.0, min_translation_threshold=0.1, minimum_points=500, agreement_threshold=0.25)#

Run pc_align against ICESat-2 and evaluate whether to keep the result.

Wraps alignment_report() with a decision step so the aligned DEM is only retained when it represents a meaningful improvement. Returns an AlignmentResult; notebook callers can inspect result.status to decide what to display.

Decision logic:

  1. If the alignment report is empty (fewer than minimum_points ICESat-2 points, or the pc_align log is missing), delete any aligned DEM file and return status="insufficient_points".

  2. Otherwise compute improvement_pct = (p50_beg - p50_end) / p50_beg * 100. If p50_end >= p50_beg or improvement_pct <= improvement_threshold_pct, delete the aligned DEM, clear self.aligned_dem_fn, and return status="no_improvement".

  3. Otherwise re-run atl06sr_to_dem_dh() so the icesat_minus_aligned_dem column is populated, and return status="success" with the aligned DEM retained.

Parameters:
  • processing_level (str, optional) – ATL06-SR processing level key to align against. Default “all”.

  • improvement_threshold_pct (float, optional) – Minimum required (p50_beg - p50_end) / p50_beg * 100 for the aligned DEM to be kept. Default 5.0.

  • min_translation_threshold (float, optional) – Forwarded to alignment_report(). Default 0.1.

  • minimum_points (int, optional) – Forwarded to alignment_report(). Default 500.

  • agreement_threshold (float, optional) – Forwarded to alignment_report(). Default 0.25.

Return type:

AlignmentResult

align_and_evaluate_planetary(max_displacement=500, improvement_threshold_pct=5.0, min_translation_threshold=0.1, minimum_points=20)#

Run pc_align against MOLA/LOLA and evaluate whether to keep it.

Mirrors align_and_evaluate() (the ICESat-2 path) but for planetary altimetry. Requires load_planetary_csv() to have been called so self.planetary_points is populated.

Decision logic:

  1. If fewer than minimum_points valid planetary points are available, return status="insufficient_points".

  2. Otherwise compute improvement_pct = (p50_beg - p50_end) / p50_beg * 100. If p50_end >= p50_beg or improvement_pct <= improvement_threshold_pct or the translation magnitude is below min_translation_threshold × DEM GSD, delete the aligned DEM and return status="no_improvement".

  3. Otherwise re-run planetary_to_dem_dh() to populate altimetry_minus_aligned_dem and return status="success".

Parameters:
  • max_displacement (float, optional) – --max-displacement for pc_align, in meters. Default 500 (ASAP-Stereo’s CTX cookbook recommendation).

  • improvement_threshold_pct (float, optional) – Minimum p50 reduction (%) required to keep the aligned DEM.

  • min_translation_threshold (float, optional) – Minimum translation magnitude as a fraction of the DEM GSD.

  • minimum_points (int, optional) – Minimum number of valid altimetry points to attempt alignment. Planetary tracks are sparser than ICESat-2, so this defaults to a much smaller number than the Earth path.

Return type:

AlignmentResult

alignment_report(processing_level='ground', minimum_points=500, agreement_threshold=0.25, write_out_aligned_dem=False, min_translation_threshold=0.1, key_for_aligned_dem=None)#

Generate alignment reports and optionally align the DEM.

Runs pc_align between the DEM and filtered ATL06-SR data for all temporal variations of a given processing level, generates reports of the alignment results, and optionally creates an aligned DEM.

Parameters:
  • processing_level (str, optional) – Base processing level to use, default is “ground”

  • minimum_points (int, optional) – Minimum number of points required for alignment, default is 500

  • agreement_threshold (float, optional) – Threshold for agreement between different temporal alignments, as a fraction of the mean shift, default is 0.25

  • write_out_aligned_dem (bool, optional) – Whether to create an aligned DEM, default is False

  • min_translation_threshold (float, optional) – Minimum translation magnitude as a fraction of the DEM GSD to warrant creating an aligned DEM, default is 0.1

  • key_for_aligned_dem (str, optional) – Which temporal filter key to use for alignment if write_out_aligned_dem is True. Default is None, which uses the processing_level value.

Returns:

Sets self.alignment_report_df and optionally self.aligned_dem_fn

Return type:

None

Notes

This method performs both the pc_align operations and analysis of the alignment results. It checks for consistency across temporal filters and only creates an aligned DEM if the translation is significant enough and consistent across temporal windows.

atl06sr_to_dem_dh(n_sigma=3)#

Calculate height differences between ATL06-SR data and DEMs.

Interpolates DEM heights at ATL06-SR point locations and calculates the difference between ICESat-2 heights and DEM heights. If an aligned DEM is available, also calculates differences against it. Outliers beyond n_sigma × NMAD from the median are removed by default.

Parameters:

n_sigma (float or None, optional) – Remove dh outliers beyond this many NMAD from the median. Default 3. Pass None to skip outlier filtering.

Returns:

Adds columns to the filtered ATL06-SR data: - dem_height: Interpolated height from the original DEM - icesat_minus_dem: Height difference (ICESat-2 - DEM) - aligned_dem_height: Interpolated height from aligned DEM (if available) - icesat_minus_aligned_dem: Height difference with aligned DEM (if available)

Return type:

None

Notes

This method performs bilinear interpolation of DEM values at ATL06-SR point locations using xarray and rioxarray. It handles coordinate system conversions automatically.

filter_esa_worldcover(filter_out='water', retain_only=None)#

Filter ATL06-SR data based on ESA WorldCover land cover classes.

Filters the data points based on ESA WorldCover classification, either by removing specific land cover types or retaining only specific types.

Parameters:
  • filter_out (str, optional) – Land cover group to filter out, default is “water”. Options: “water”, “snow_ice”, “trees”, “low_vegetation”, “built_up”

  • retain_only (str or None, optional) – If specified, retain only points matching this land cover group, default is None. Same options as filter_out.

Returns:

Results are stored in the class attributes

Return type:

None

Notes

This method uses the ESA WorldCover land cover classification (see WORLDCOVER_NAMES), which was sampled when requesting the ATL06-SR data.

filter_outliers(column='icesat_minus_dem', n_sigma=3)#

Remove dh outliers beyond n_sigma × standard deviation from the mean.

Parameters:
  • column (str, optional) – Column to filter on, default "icesat_minus_dem".

  • n_sigma (float, optional) – Number of standard deviations to allow, default 3.

generic_temporal_filter_atl06sr(select_years=None, select_months=None, select_days=None)#

Apply custom temporal filters to ATL06-SR data.

Filters the data based on specific years, months, or days, allowing for custom temporal filtering not covered by the predefined filters.

Parameters:
  • select_years (list or None, optional) – Years to keep (e.g., [2019, 2020]), default is None

  • select_months (list or None, optional) – Months to keep (1-12), default is None

  • select_days (list or None, optional) – Days of month to keep (1-31), default is None

Returns:

Results are filtered in-place in the class attributes

Return type:

None

Notes

This method modifies the existing filtered data rather than creating new entries with different keys. At least one of the filter parameters should be provided for the method to have any effect.

histogram(key='all', title='Histogram', plot_aligned=False, save_dir=None, fig_fn=None)#

Plot histograms of height differences between ATL06-SR data and DEMs.

Creates histograms of the height differences between ICESat-2 ATL06-SR data and DEMs, with statistics including median and normalized median absolute deviation (NMAD).

Parameters:
  • key (str, optional) – Processing level to plot, default is “all”

  • title (str, optional) – Plot title, default is “Histogram”

  • plot_aligned (bool, optional) – Whether to include differences with aligned DEM, default is False

  • save_dir (str or None, optional) – Directory to save figure, default is None (don’t save)

  • fig_fn (str or None, optional) – Filename for saved figure, default is None

Returns:

Displays the plot and optionally saves it

Return type:

None

Notes

If the height differences haven’t been calculated yet, this method calls atl06sr_to_dem_dh() to calculate them. NMAD is a robust measure of dispersion that is less sensitive to outliers than standard deviation, calculated as 1.4826 * median(abs(x - median(x))).

histogram_by_landcover(key='all', top_n=4, title='ICESat-2 ATL06-SR vs DEM', xlim=None, plot_aligned=False, save_dir=None, fig_fn=None)#

Plot histogram of dh with per-landcover-class statistics.

Creates a histogram of the height differences between ICESat-2 ATL06-SR data and the DEM, with a text annotation showing overall and per-landcover-class statistics (count, median, NMAD).

When plot_aligned=True and an aligned DEM is available, overlays the pre- and post-alignment distributions and renders two vertically stacked stats text boxes whose outline colors match the bar colors (color serves as the legend).

Parameters:
  • key (str, optional) – Processing level key, default is “all”

  • top_n (int, optional) – Number of top landcover classes to report, default is 4

  • title (str, optional) – Plot title, default is “ICESat-2 ATL06-SR vs DEM”

  • xlim (tuple or None, optional) – Symmetric x-axis limits as (min, max). If None, uses ±3σ range (data is already 3σ-filtered in atl06sr_to_dem_dh).

  • plot_aligned (bool, optional) – Whether to overlay the aligned-DEM distribution alongside the unaligned one. Requires self.aligned_dem_fn and the icesat_minus_aligned_dem column. Default is False.

  • save_dir (str or None, optional) – Directory to save figure, default is None

  • fig_fn (str or None, optional) – Filename for saved figure, default is None

histogram_planetary_to_dem(save_dir=None, fig_fn=None, title=None, plot_aligned=False)#

Histogram of planetary altimetry vs DEM height differences.

Parameters:
  • save_dir (str or None, optional) – Directory to save figure.

  • fig_fn (str or None, optional) – Filename for saved figure.

  • title (str or None, optional) – Custom plot title. Auto-detected if None.

  • plot_aligned (bool, optional) – Overlay the post-alignment dh distribution on the same axes. Requires that pc_align has been run successfully.

load_planetary_csv(csv_path)#

Load LOLA or MOLA altimetry data from a GDS topo CSV file.

The CSV is obtained via the request_planetary_altimetry CLI tool, which submits an async query to the ODE GDS API and emails the user a download link. The user downloads and unzips the result, then passes the *_topo_csv.csv file here.

Automatically selects the LOLA or MOLA parser based on the DEM’s planetary body.

Parameters:

csv_path (str) – Path to a *_topo_csv.csv file from the ODE GDS.

mapview_plot_atl06sr_to_dem(key='all', clim=None, plot_aligned=False, save_dir=None, fig_fn=None, map_crs=None, **ctx_kwargs)#

Plot height differences between ATL06-SR data and DEMs.

Creates a map visualization of the height differences between ICESat-2 ATL06-SR data and either the original or aligned DEM.

Parameters:
  • key (str, optional) – Processing level to plot, default is “all”

  • clim (tuple or None, optional) – Color limits as (min, max), default is None (auto)

  • plot_aligned (bool, optional) – Whether to plot differences with aligned DEM, default is False

  • save_dir (str or None, optional) – Directory to save figure, default is None (don’t save)

  • fig_fn (str or None, optional) – Filename for saved figure, default is None

  • map_crs (str or None, optional) – Coordinate reference system for mapping, default is None (use DEM’s CRS)

  • **ctx_kwargs (dict, optional) – Additional arguments for contextily basemap

Returns:

Displays the plot and optionally saves it

Return type:

None

Notes

If the height differences haven’t been calculated yet, this method calls atl06sr_to_dem_dh() to calculate them. The plot uses a divergent colormap (RdBu) to highlight positive and negative differences.

mapview_plot_planetary_to_dem(clim=None, save_dir=None, fig_fn=None, title=None, plot_aligned=False)#

Map view of planetary altimetry vs DEM height differences.

Plots the DEM hillshade as background with altimetry dh points overlaid using a divergent colourmap. When plot_aligned=True and self.aligned_dem_fn is set, renders pre/post panels side by side.

Parameters:
  • clim (tuple or None, optional) – Colour limits (min, max) for dh. Default auto (symmetric ±|max| around zero).

  • save_dir (str or None, optional) – Directory to save figure.

  • fig_fn (str or None, optional) – Filename for saved figure.

  • title (str or None, optional) – Custom plot title. Auto-detected if None.

  • plot_aligned (bool, optional) – Add a second panel showing dh against the aligned DEM. Requires that pc_align has been run successfully.

planetary_to_dem_dh(n_sigma=3)#

Compute height differences between planetary altimetry and DEM.

Reprojects self.planetary_points to the DEM CRS, interpolates DEM heights at altimetry locations, and computes the difference altimetry_minus_dem = height - dem_height. When self.aligned_dem_fn is set, also populates aligned_dem_height and altimetry_minus_aligned_dem so pre/post-alignment plots can share a single sample. Outliers beyond n_sigma × std from the mean (computed on the unaligned dh) are removed by default.

Parameters:
  • n_sigma (float or None, optional) – Remove dh outliers beyond this many standard deviations from the mean. Default 3. Pass None to skip outlier filtering.

  • self.planetary_points. (The results are stored as new columns on)

plot_atl06sr(key='all', plot_beams=False, plot_dem=False, column_name='h_mean', cbar_label='Height above datum (m)', title='ICESat-2 ATL06-SR', clim=None, symm_clim=False, cmap='inferno', map_crs='EPSG:4326', figsize=(6, 4), save_dir=None, fig_fn=None, **ctx_kwargs)#

Plot ATL06-SR data on a map with customizable options.

Creates a map view of ATL06-SR data with options to color by various attributes, highlight different laser beams, overlay on the DEM, and add contextual basemaps.

Parameters:
  • key (str, optional) – Processing level to plot, default is “all”

  • plot_beams (bool, optional) – Whether to color points by ICESat-2 beam, default is False

  • plot_dem (bool, optional) – Whether to plot the DEM as a background, default is False

  • column_name (str, optional) – Column to use for point coloring, default is “h_mean”

  • cbar_label (str, optional) – Colorbar label, default is “Height above datum (m)”

  • title (str, optional) – Plot title, default is “ICESat-2 ATL06-SR”

  • clim (tuple or None, optional) – Color limits as (min, max), default is None (auto)

  • symm_clim (bool, optional) – Whether to use symmetric color limits, default is False

  • cmap (str, optional) – Matplotlib colormap, default is “inferno”

  • map_crs (str, optional) – Coordinate reference system for mapping, default is “EPSG:4326”

  • figsize (tuple, optional) – Figure size as (width, height), default is (6, 4)

  • save_dir (str or None, optional) – Directory to save figure, default is None (don’t save)

  • fig_fn (str or None, optional) – Filename for saved figure, default is None

  • **ctx_kwargs (dict, optional) – Additional arguments for contextily basemap

Returns:

Displays the plot and optionally saves it

Return type:

None

Notes

When plot_beams is True, points are colored by ICESat-2 laser spot number, with strong beams (1, 3, 5) in darker colors and weak beams (2, 4, 6) in lighter colors.

plot_atl06sr_dem_profile(key='all', rgt=None, cycle=None, spot=None, plot_aligned=False, save_dir=None, fig_fn=None)#

Plot elevation profile comparing ICESat-2 and DEM along the best track.

Creates a 2×2 figure with the profile stack on the left and a map view spanning the full height on the right: - Top-left: Absolute elevation profile (DEM, COP30, ICESat-2) - Bottom-left: Height difference profile (ICESat-2 minus DEM)

(shares x-axis with top-left, no vertical space between them)

  • Right column: DEM hillshade map with the full track and segment extents, spanning the full vertical height

Parameters:
  • key (str, optional) – Processing level key, default is “all”

  • rgt (int or None, optional) – Reference ground track (auto-selected if None)

  • cycle (int or None, optional) – Cycle number (auto-selected if None)

  • spot (int or None, optional) – Spot number (auto-selected if None)

  • plot_aligned (bool, optional) – Whether to also plot the aligned DEM profile, default is False

  • save_dir (str or None, optional) – Directory to save figure, default is None

  • fig_fn (str or None, optional) – Filename for saved figure, default is None

plot_atl06sr_time_stamps(key='all', title='ICESat-2 ATL06-SR Time Stamps', cmap='inferno', map_crs='EPSG:4326', figsize=(15, 10), save_dir=None, fig_fn=None, **ctx_kwargs)#

Plot ATL06-SR data for different temporal filters.

Creates a 2x2 grid of plots showing ATL06-SR data for different temporal filters (unfiltered, 15-day, 45-day, and seasonal) colored by height.

Parameters:
  • key (str, optional) – Base processing level to plot, default is “all”

  • title (str, optional) – Plot title, default is “ICESat-2 ATL06-SR Time Stamps”

  • cmap (str, optional) – Matplotlib colormap for elevation, default is “inferno”

  • map_crs (str, optional) – Coordinate reference system for mapping, default is “EPSG:4326”

  • figsize (tuple, optional) – Figure size as (width, height), default is (15, 10)

  • save_dir (str or None, optional) – Directory to save figure, default is None (don’t save)

  • fig_fn (str or None, optional) – Filename for saved figure, default is None

  • **ctx_kwargs (dict, optional) – Additional arguments for contextily basemap

Returns:

Displays the plot and optionally saves it

Return type:

None

Notes

This method requires the filtered data to have been created using the predefined_temporal_filter_atl06sr method for the temporal variations to be available.

plot_best_worst_segments(key='all', rgt=None, cycle=None, spot=None, plot_aligned=False, save_dir=None, fig_fn=None)#

Plot 1 km segments with better and worse agreement as a 1×2 figure.

Creates a single-row, 2-column figure: - Column 1: Better agreement segment (lowest score) - Column 2: Worse agreement segment (highest score)

Segment score is 3·|median(dh)| + NMAD(dh) (see _find_best_worst_segments). Segment selection is based on the unaligned dh so best/worst segments remain comparable across the pre- and post-alignment variants of this plot.

Parameters:
  • key (str, optional) – Processing level key, default is “all”

  • rgt (int or None, optional) – Reference ground track (auto-selected if None)

  • cycle (int or None, optional) – Cycle number (auto-selected if None)

  • spot (int or None, optional) – Spot number (auto-selected if None)

  • plot_aligned (bool, optional) – Whether to overlay the aligned DEM heights and include aligned Median/NMAD in each segment title. Requires self.aligned_dem_fn and the aligned_dem_height / icesat_minus_aligned_dem columns. Default False.

  • save_dir (str or None, optional) – Directory to save figure, default is None

  • fig_fn (str or None, optional) – Filename for saved figure, default is None

predefined_temporal_filter_atl06sr(date=None)#

Apply predefined temporal filters to ATL06-SR data.

Creates multiple temporally filtered versions of the ATL06-SR data based on a reference date, including: - 15-day window around the date - 45-day window around the date - 91-day window around the date - Seasonal filter (same season as the reference date)

Parameters:

date (datetime or None, optional) – Reference date for filtering. If None, extracts date from stereopair metadata, default is None

Returns:

Results are stored in the class attributes with new keys indicating the temporal filter applied

Return type:

None

Notes

This method is particularly useful for DEM validation and alignment, as it provides multiple temporal windows to analyze the stability of the terrain and to identify optimal temporal windows for alignment.

request_atl06sr_multi_processing(processing_levels=['all', 'ground', 'canopy', 'top_of_canopy'], res=20, len=40, ats=20, cnt=10, maxi=6, h_sigma_quantile=1.0, save_to_parquet=False, filename='atl06sr', region=None, time_range='all', scene_date=None, time_buffer_days=365, t0=None, t1=None)#

Request ICESat-2 ATL06-SR data for multiple processing levels.

Downloads ATL06-SR data from the SlideRule API for specified processing levels (surface types), with options to filter and save the results. Each processing level targets different surface types like ground, canopy, etc.

Parameters:
  • processing_levels (list, optional) – List of processing levels to request, default is [“all”, “ground”, “canopy”, “top_of_canopy”]

  • res (int, optional) – ATL06-SR segment resolution in meters, default is 20

  • len (int, optional) – ATL06-SR segment length in meters, default is 40

  • ats (int, optional) – Along-track sigma, default is 20

  • cnt (int, optional) – Minimum number of photons for segment, default is 10

  • maxi (int, optional) – Maximum iterations for surface fit, default is 6

  • h_sigma_quantile (float, optional) – Quantile for filtering by h_sigma, default is 1.0

  • save_to_parquet (bool, optional) – Whether to save results to parquet files, default is False

  • filename (str, optional) – Base filename for saved data, default is “atl06sr”

  • region (list or None, optional) – Region bounds as [minx, miny, maxx, maxy] in lat/lon, default is None (derived from DEM)

  • time_range (str, optional) – "all" (default) requests all ICESat-2 data from mission start to present. "buffered" activates time filtering via the cascade: t0/t1 > scene_date ± time_buffer_days > XML metadata ± time_buffer_days > fall back to all.

  • scene_date (str or datetime-like, optional) – Scene acquisition date, used when time_range="buffered". If None, auto-detected from stereopair XML metadata.

  • time_buffer_days (int, optional) – Days before/after scene_date defining the time window, default is 365

  • t0 (str or datetime-like, optional) – Explicit start date (e.g. “2020-01-01”), used when time_range="buffered". Overrides scene_date.

  • t1 (str or datetime-like, optional) – Explicit end date (e.g. “2024-12-31”). Defaults to present if only t0 is provided.

Returns:

Results are stored in the class attributes

Return type:

None

Notes

This method makes SlideRule API calls which require internet connectivity. It also samples ESA WorldCover data for land cover classification of the ICESat-2 data points. The method includes filtering to improve data quality by removing points with high uncertainty and from early mission cycles.

sample_esa_worldcover()#

Sample ESA WorldCover 10m values into every filtered processing level that doesn’t already have them.

Normally not needed when data is requested via request_atl06sr_multi_processing, which samples WorldCover automatically before caching. Use this when working with manually-loaded data.

to_csv_for_pc_align(key='ground', filename_prefix='atl06sr_for_pc_align')#

Export ATL06-SR data to CSV format for use with pc_align.

Creates a CSV file from the filtered ATL06-SR data that is formatted for use with the ASP pc_align tool, containing longitude, latitude, and height above datum.

Parameters:
  • key (str, optional) – Processing level to export, default is “ground”

  • filename_prefix (str, optional) – Prefix for the output CSV file, default is “atl06sr_for_pc_align”

Returns:

Path to the created CSV file

Return type:

str

Notes

The ASP pc_align tool requires input data in a specific format. This method converts the ATL06-SR GeoDataFrame to the required CSV format with columns for longitude, latitude, and height.

to_csv_for_pc_align_planetary(filename_prefix='planetary_for_pc_align')#

Export self.planetary_points to a CSV for pc_align.

Writes columns lon, lat, radius_m (planetary radius from the body center, in meters). Used as the planetary_csv argument to asp_plot.alignment.Alignment.pc_align_dem_to_planetary_csv().

Parameters:

filename_prefix (str, optional) – Prefix for the output CSV filename. Saved in self.directory.

Returns:

Absolute path to the created CSV file.

Return type:

str

aligned_dem_fn = None#
atl06sr_processing_levels#
atl06sr_processing_levels_filtered#
dem_fn#
directory#
planetary_points = None#
asp_plot.altimetry.gds_query_async(query_type, bounds, results_code, email=None, **extra_params)#

Submit an async query to the ODE GDS REST API.

Parameters:
  • query_type (str) – GDS query type, e.g. "lolardr" or "molapedr".

  • bounds (dict) – Dictionary with westernlon, easternlon, minlat, maxlat keys.

  • results_code (str) – GDS results format code (e.g. "u" for LOLA, "v" for MOLA).

  • email (str or None, optional) – Email for notification when query finishes.

  • **extra_params – Additional GDS query parameters (e.g. channel="ttttt").

Returns:

Job ID for polling.

Return type:

str

asp_plot.altimetry.GDS_BASE_URL = 'https://oderest.rsl.wustl.edu/livegds'#
asp_plot.altimetry.ICESAT2_MISSION_START#
asp_plot.altimetry.MARS_IAU_SPHERE_RADIUS = 3396190.0#
asp_plot.altimetry.MOON_IAU_SPHERE_RADIUS = 1737400.0#
asp_plot.altimetry.WORLDCOVER_NAMES#
asp_plot.altimetry.logger#