API Reference

This page lists all available functions and classes in OGGM. It is a hard work to keep everything up-to-date, so don’t hesitate to let us know (see Get in touch) if something’s missing, or help us (see Contributing to OGGM) to write a better documentation!

Workflow

Tools to set-up and run OGGM.

cfg.initialize

Read the configuration file containing the run’s parameters.

cfg.set_logging_config

Set the global logger parameters.

cfg.set_intersects_db

Set the glacier intersection database for OGGM to use.

cfg.reset_working_dir

Deletes the content of the working directory.

workflow.init_glacier_regions

Initializes the list of Glacier Directories for this run.

workflow.execute_entity_task

Execute a task on gdirs.

workflow.gis_prepro_tasks

Shortcut function: run all flowline preprocessing tasks.

workflow.climate_tasks

Shortcut function: run all climate related tasks.

workflow.inversion_tasks

Shortcut function: run all ice thickness inversion tasks.

workflow.merge_glacier_tasks

Shortcut function: run all tasks to merge tributaries to a main glacier

utils.compile_glacier_statistics

Gather as much statistics as possible about a list of glaciers.

utils.compile_run_output

Merge the output of the model runs of several gdirs into one file.

utils.compile_climate_input

Merge the climate input files in the glacier directories into one file.

utils.compile_task_log

Gathers the log output for the selected task(s)

utils.copy_to_basedir

Copies the glacier directories and their content to a new location.

utils.gdir_to_tar

Writes the content of a glacier directory to a tar file.

utils.base_dir_to_tar

Merge the directories into 1000 bundles as tar files.

Troubleshooting

utils.show_versions

Prints the OGGM version and other system information.

Input/Output

utils.get_demo_file

Returns the path to the desired OGGM-sample-file.

utils.get_rgi_dir

Path to the RGI directory.

utils.get_rgi_region_file

Path to the RGI region file.

utils.get_rgi_glacier_entities

Get a list of glacier outlines selected from their RGI IDs.

utils.get_rgi_intersects_dir

Path to the RGI directory containing the intersect files.

utils.get_rgi_intersects_region_file

Path to the RGI regional intersect file.

utils.get_rgi_intersects_entities

Get a list of glacier intersects selected from their RGI IDs.

utils.get_cru_file

Returns a path to the desired CRU baseline climate file.

utils.get_histalp_file

Returns a path to the desired HISTALP baseline climate file.

utils.get_ref_mb_glaciers

Get the list of glaciers we have valid mass balance measurements for.

Entity tasks

Entity tasks are tasks which are applied on single glaciers individually and do not require information from other glaciers (this encompasses the majority of OGGM’s tasks). They are parallelizable.

tasks.define_glacier_region

Very first task: define the glacier’s local grid.

tasks.glacier_masks

Makes a gridded mask of the glacier outlines and topography.

tasks.compute_centerlines

Compute the centerlines following Kienholz et al., (2014).

tasks.initialize_flowlines

Computes more physical Inversion Flowlines from geometrical Centerlines

tasks.compute_downstream_line

Computes the Flowline along the unglaciated downstream topography

tasks.compute_downstream_bedshape

The bedshape obtained by fitting a parabola to the line’s normals.

tasks.catchment_area

Compute the catchment areas of each tributary line.

tasks.catchment_intersections

Computes the intersections between the catchments.

tasks.catchment_width_geom

Compute geometrical catchment widths for each point of the flowlines.

tasks.catchment_width_correction

Corrects for NaNs and inconsistencies in the geometrical widths.

tasks.process_cru_data

Processes and writes the CRU baseline climate data for this glacier.

tasks.process_histalp_data

Processes and writes the HISTALP baseline climate data for this glacier.

tasks.process_custom_climate_data

Processes and writes the climate data from a user-defined climate file.

tasks.process_gcm_data

Applies the anomaly method to GCM climate data

tasks.process_cesm_data

Processes and writes CESM climate data for this glacier.

tasks.process_cmip5_data

Read, process and store the CMIP5 climate data data for this glacier.

tasks.local_t_star

Compute the local t* and associated glacier-wide mu*.

tasks.mu_star_calibration

Compute the flowlines’ mu* and the associated apparent mass-balance.

tasks.apparent_mb_from_linear_mb

Compute apparent mb from a linear mass-balance assumption (for testing).

tasks.glacier_mu_candidates

Computes the mu candidates, glacier wide.

tasks.prepare_for_inversion

Prepares the data needed for the inversion.

tasks.mass_conservation_inversion

Compute the glacier thickness along the flowlines

tasks.filter_inversion_output

Filters the last few grid point whilst conserving total volume.

tasks.distribute_thickness_per_altitude

Compute a thickness map by redistributing mass along altitudinal bands.

tasks.distribute_thickness_interp

Compute a thickness map by interpolating between centerlines and border.

tasks.init_present_time_glacier

Merges data from preprocessing tasks.

tasks.run_random_climate

Runs the random mass-balance model for a given number of years.

tasks.run_constant_climate

Runs the constant mass-balance model for a given number of years.

tasks.run_from_climate_data

Runs a glacier with climate input from e.g.

Global tasks

Global tasks are tasks which are run on a set of glaciers (most often: all glaciers in the current run). They are not parallelizable at the user level but might use multiprocessing internally.

tasks.compute_ref_t_stars

Detects the best t* for the reference glaciers and writes them to disk

tasks.compile_glacier_statistics

Gather as much statistics as possible about a list of glaciers.

tasks.compile_run_output

Merge the output of the model runs of several gdirs into one file.

tasks.compile_climate_input

Merge the climate input files in the glacier directories into one file.

tasks.compile_task_log

Gathers the log output for the selected task(s)

Classes

Listed here are the classes which are relevant at the API level (i.e. classes which are used and re-used across modules and tasks).

GlacierDirectory

Organizes read and write access to the glacier’s files.

Centerline

A Centerline has geometrical and flow rooting properties.

Flowline

The is the input flowline for the model.

Glacier directories

Glacier directories (see also: GlacierDirectory) are folders on disk which store the input and output data for a single glacier during an OGGM run. The data are on disk to be persistent, i.e. they won’t be deleted unless you ask OGGM to. You can start a run from an existing directory, avoiding to re-do unnecessary computations.

Initialising a glacier directory

The easiest way to initialize a glacier directory is to start from a pre-processed state available on the OGGM servers:

In [1]: import os

In [2]: import oggm

In [3]: from oggm import cfg, workflow

In [4]: from oggm.utils import gettempdir

In [5]: cfg.initialize()  # always initialize before an OGGM task

# The working directory is where OGGM will store the run's data
In [6]: cfg.PATHS['working_dir'] = os.path.join(gettempdir(), 'Docs_GlacierDir')

In [7]: gdirs = workflow.init_glacier_regions('RGI60-11.00897', from_prepro_level=1,
   ...:                                       prepro_border=10)
   ...: 

In [8]: gdir = gdirs[0]  # init_glacier_regions always returns a list

We just downloaded the minimal input for a glacier directory. The GlacierDirectory object contains the glacier metadata:

In [9]: gdir
Out[9]: 
<oggm.GlacierDirectory>
  RGI id: RGI60-11.00897
  Region: 11: Central Europe
  Subregion: 11-01: Alps                            
  Name: Hintereisferner
  Glacier type: Glacier
  Terminus type: Land-terminating
  Area: 8.036 km2
  Lon, Lat: (10.7584, 46.8003)
  Grid (nx, ny): (139, 94)
  Grid (dx, dy): (50.0, -50.0)

In [10]: gdir.rgi_id
Out[10]: 'RGI60-11.00897'

It also points to a specific directory in disk, where the files are written to:

In [11]: gdir.dir
Out[11]: '/tmp/OGGM/Docs_GlacierDir/per_glacier/RGI60-11/RGI60-11.00/RGI60-11.00897'

In [12]: os.listdir(gdir.dir)
Out[12]: 
['dem.tif',
 'diagnostics.json',
 'dem_source.txt',
 'glacier_grid.json',
 'outlines.tar.gz.properties',
 'outlines.tar.gz',
 'intersects.tar.gz',
 'log.txt']

Users usually don’t have to care about where the data is located. GlacierDirectory objects help you to get this information:

In [13]: fdem = gdir.get_filepath('dem')

In [14]: fdem
Out[14]: '/tmp/OGGM/Docs_GlacierDir/per_glacier/RGI60-11/RGI60-11.00/RGI60-11.00897/dem.tif'

In [15]: import xarray as xr

In [16]: xr.open_rasterio(fdem).plot(cmap='terrain')
Out[16]: <matplotlib.collections.QuadMesh at 0x7f74d2a70588>
_images/plot_gdir_dem.png

This persistence on disk allows for example to continue a workflow that has been previously interrupted. Initialising a GlacierDirectory from a non-empty folder won’t erase its content:

In [17]: gdir = oggm.GlacierDirectory('RGI60-11.00897')

In [18]: os.listdir(gdir.dir)  # the directory still contains the data
Out[18]: 
['dem.tif',
 'diagnostics.json',
 'dem_source.txt',
 'glacier_grid.json',
 'outlines.tar.gz.properties',
 'outlines.tar.gz',
 'intersects.tar.gz',
 'log.txt']

cfg.BASENAMES

This is a list of the files that can be found in the glacier directory or its divides. These data files and their names are standardized and listed in the oggm.cfg module. If you want to implement your own task you’ll have to add an entry to this file too.

calving_loop.csv

A csv file containing the output of each iteration of the find_inversion_calving task loop.

calving_output.pkl

Calving output (deprecated).

catchments_intersects.shp

The intersections between cathments (shapefile) in the local map projection (Transverse Mercator).

centerlines.pkl

A list of oggm.Centerline instances, sorted by flow order.

climate_info.pkl

Some information (dictionary) about the climate data and the mass balance parameters for this glacier.

climate_monthly.nc

The monthly climate timeseries stored in a netCDF file.

dem.tif

A geotiff file containing the DEM (reprojected into the local grid).This DEM is not smoothed or gap filles, and is the closest to the original DEM source.

dem_source.txt

A text file with the source of the topo file (GIMP, SRTM, …).

diagnostics.json

A dictionary containing runtime diagnostics useful for debugging.

downstream_line.pkl

A dictionary containing the downsteam line geometry as well as the bed shape computed from a parabolic fit.

flowline_catchments.shp

Each flowline has a catchment area computed from flow routing algorithms: this shapefile stores the catchment outlines (in the local map projection (Transverse Mercator).

gcm_data.nc

The monthly GCM climate timeseries stored in a netCDF file.

geometries.pkl

A dictionary containing the shapely.Polygons of a glacier. The “polygon_hr” entry contains the geometry transformed to the local grid in (i, j) coordinates, while the “polygon_pix” entry contains the geometries transformed into the coarse grid (the i, j elements are integers). The “polygon_area” entry contains the area of the polygon as computed by Shapely. The “catchment_indices” entrycontains a list of len n_centerlines, each element containing a numpy array of the indices in the glacier grid which represent the centerlines catchment area.

glacier_grid.json

A salem.Grid handling the georeferencing of the local grid.

gridded_data.nc

A netcdf file containing several gridded data variables such as topography, the glacier masks, the interpolated 2D glacier bed, and more.

hypsometry.csv

A hypsometry file computed by OGGM and provided in the same format as the RGI (useful for diagnostics).

intersects.shp

The glacier intersects in the local map projection (Transverse Mercator).

inversion_flowlines.pkl

A “better” version of the Centerlines, now on a regular spacing i.e., not on the gridded (i, j) indices. The tails of the tributaries are cut out to make more realistic junctions. They are now “1.5D” i.e., with a width.

inversion_input.pkl

List of dicts containing the data needed for the inversion.

inversion_output.pkl

List of dicts containing the output data from the inversion.

inversion_params.pkl

Dict of fs and fd as computed by the inversion optimisation.

linear_mb_params.pkl

When using a linear mass-balance for the inversion, this dict stores the optimal ela_h and grad.

local_mustar.json

A dict containing the glacier’s t*, bias, and the flowlines’ mu*

model_diagnostics.nc

A netcdf file containing the model diagnostics (volume, mass-balance, length…).

model_flowlines.pkl

List of flowlines ready to be run by the model.

model_run.nc

A netcdf file containing enough information to reconstruct the entire flowline glacier along the run (can be data expensive).

outlines.shp

The glacier outlines in the local map projection (Transverse Mercator).

Mass-balance

Mass-balance models found in the core.massbalance module.

They follow the MassBalanceModel() interface. Here is a quick summary of the units and conventions used by all models:

Units

The computed mass-balance is in units of [m ice s-1] (“meters of ice per second”), unless otherwise specified (e.g. for the utility function get_specific_mb). The conversion from the climatic mass-balance ([kg m-2 s-1] ) therefore assumes an ice density given by cfg.PARAMS['ice_density'] (currently: 900 kg m-3).

Time

The time system used by OGGM is a simple “fraction of year” system, where the floating year can be used for conversion to months and years:

In [19]: from oggm.utils import floatyear_to_date, date_to_floatyear

In [20]: date_to_floatyear(1982, 12)
Out[20]: 1982.9166666666667

In [21]: floatyear_to_date(1.2)
Out[21]: (1, 3)

Interface

MassBalanceModel

Common logic for the mass balance models.

MassBalanceModel.get_monthly_mb

Monthly mass-balance at given altitude(s) for a moment in time.

MassBalanceModel.get_annual_mb

Like self.get_monthly_mb(), but for annual MB.

MassBalanceModel.get_specific_mb

Specific mb for this year and a specific glacier geometry.

MassBalanceModel.get_ela

Compute the equilibrium line altitude for this year

Models

LinearMassBalance

Constant mass-balance as a linear function of altitude.

PastMassBalance

Mass balance during the climate data period.

ConstantMassBalance

Constant mass-balance during a chosen period.

RandomMassBalance

Random shuffle of all MB years within a given time period.

UncertainMassBalance

Adding uncertainty to a mass balance model.

MultipleFlowlineMassBalance

Handle mass-balance at the glacier level instead of flowline level.