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.initialize_minimal Same as initialise() but without requiring any download of data.
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.
cfg.add_to_basenames Add an entry to the list of BASENAMES.
workflow.init_glacier_directories Initializes the list of Glacier Directories for this run.
workflow.execute_entity_task Execute a task on gdirs.
workflow.gis_prepro_tasks Run all flowline preprocessing tasks on a list of glaciers.
workflow.download_ref_tstars Downloads and copies the reference list of t* to the working directory.
workflow.climate_tasks Run all climate related entity tasks on a list of glaciers.
workflow.inversion_tasks Run all ice thickness inversion tasks on a list of glaciers.
workflow.merge_glacier_tasks Shortcut function: run all tasks to merge tributaries to a main glacier
workflow.calibrate_inversion_from_consensus Fit the total volume of the glaciers to the 2019 consensus estimate.
workflow.match_regional_geodetic_mb Regional shift of the mass-balance residual to match observations.

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_ref_mb_glaciers Get the list of glaciers we have valid mass balance measurements for.
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.
global_tasks.write_centerlines_to_shape Write the centerlines in a shapefile.
global_tasks.compile_glacier_statistics Gather as much statistics as possible about a list of glaciers.
global_tasks.compile_run_output Merge the output of the model runs of several gdirs into one file.
global_tasks.compile_climate_input Merge the climate input files in the glacier directories into one file.
global_tasks.compile_task_log Gathers the log output for the selected task(s)
global_tasks.compile_task_time Gathers the time needed for the selected task(s) to run
global_tasks.compile_fixed_geometry_mass_balance Compiles a table of specific mass-balance timeseries for all glaciers.
global_tasks.compile_climate_statistics Gather as much statistics as possible about a list of glaciers.

OGGM Shop

shop.cru.get_cru_file Returns a path to the desired CRU baseline climate file.
shop.cru.process_cru_data Processes and writes the CRU baseline climate data for this glacier.
shop.cru.process_dummy_cru_file Create a simple baseline climate file for this glacier - for testing!
shop.ecmwf.get_ecmwf_file Returns a path to the desired ECMWF baseline climate file.
shop.ecmwf.process_ecmwf_data Processes and writes the ECMWF baseline climate data for this glacier.
shop.histalp.get_histalp_file Returns a path to the desired HISTALP baseline climate file.
shop.histalp.process_histalp_data Processes and writes the HISTALP baseline climate data for this glacier.
shop.gcm_climate.process_gcm_data Applies the anomaly method to GCM climate data
shop.gcm_climate.process_cesm_data Processes and writes CESM climate data for this glacier.
shop.gcm_climate.process_cmip_data Read, process and store the CMIP5 and CMIP6 climate data for this glacier.
shop.bedtopo.add_consensus_thickness Add the consensus thickness estimate to the gridded_data file.
shop.its_live.velocity_to_gdir Reproject the its_live files to the given glacier directory.
shop.rgitopo.init_glacier_directories_from_rgitopo Initialize a glacier directory from an RGI-TOPO DEM.
shop.rgitopo.select_dem_from_dir Select a DEM from the ones available in an RGI-TOPO glacier directory.
shop.rgitopo.dem_quality_check Run a simple quality check on the rgitopo DEMs

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 after initialization: define the glacier’s local grid.
tasks.process_dem Reads the DEM from the tiff, attempts to fill voids and apply smooth.
tasks.glacier_masks Makes a gridded mask of the glacier outlines that can be used by OGGM.
tasks.simple_glacier_masks Compute glacier masks based on much simpler rules than OGGM’s default.
tasks.rasterio_glacier_mask Writes a 1-0 glacier mask GeoTiff with the same dimensions as dem.tif
tasks.gridded_attributes Adds attributes to the gridded file, useful for thickness interpolation.
tasks.gridded_mb_attributes Adds mass-balance related attributes to the gridded data file.
tasks.gridded_data_var_to_geotiff Writes a NetCDF variable to a georeferenced geotiff file.
tasks.compute_centerlines Compute the centerlines following Kienholz et al., (2014).
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.initialize_flowlines Computes more physical Inversion Flowlines from geometrical Centerlines
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.terminus_width_correction Sets a new value for the terminus width.
tasks.elevation_band_flowline Compute “squeezed” or “collapsed” glacier flowlines from Huss 2012.
tasks.fixed_dx_elevation_band_flowline Converts the “collapsed” flowline into a regular “inversion flowline”.
tasks.glacier_mu_candidates Computes the mu candidates, glacier wide.
tasks.process_climate_data Adds the selected climate data to this glacier directory.
tasks.process_custom_climate_data Processes and writes the climate data from a user-defined climate file.
tasks.historical_delta_method Applies the anomaly method to historical climate data.
tasks.historical_climate_qc Check the “quality” of the baseline climate data and correct if needed.
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.apparent_mb_from_any_mb Compute apparent mb from an arbitrary mass-balance profile.
tasks.fixed_geometry_mass_balance Computes the mass-balance with climate input from e.g.
tasks.process_cru_data Processes and writes the CRU baseline climate data for this glacier.
tasks.process_dummy_cru_file Create a simple baseline climate file for this glacier - for testing!
tasks.process_histalp_data Processes and writes the HISTALP baseline climate data for this glacier.
tasks.process_ecmwf_data Processes and writes the ECMWF baseline climate data for this glacier.
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 Renamed to process_cmip_data.
tasks.process_cmip_data Read, process and store the CMIP5 and CMIP6 climate data for this glacier.
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 points after the physically-based inversion.
tasks.get_inversion_volume Small utility task to get to the volume od all glaciers.
tasks.compute_velocities Surface velocities along the flowlines from inverted ice thickness.
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.find_inversion_calving Optimized search for a calving flux compatible with the bed inversion.
tasks.find_inversion_calving_from_any_mb Optimized search for a calving flux compatible with the bed inversion.
tasks.init_present_time_glacier Merges data from preprocessing tasks.
tasks.flowline_model_run Runs a model simulation with the default time stepping scheme.
tasks.run_random_climate Runs the random mass-balance model for a given number of years.
tasks.run_from_climate_data Runs a glacier with climate input from e.g.
tasks.run_constant_climate Runs the constant mass-balance model for a given number of years.
tasks.copy_to_basedir Copies the glacier directories and their content to a new location.
tasks.gdir_to_tar Writes the content of a glacier directory to a tar file.

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.

global_tasks.gis_prepro_tasks Run all flowline preprocessing tasks on a list of glaciers.
global_tasks.climate_tasks Run all climate related entity tasks on a list of glaciers.
global_tasks.inversion_tasks Run all ice thickness inversion tasks on a list of glaciers.
global_tasks.calibrate_inversion_from_consensus Fit the total volume of the glaciers to the 2019 consensus estimate.
global_tasks.match_regional_geodetic_mb Regional shift of the mass-balance residual to match observations.
global_tasks.merge_glacier_tasks Shortcut function: run all tasks to merge tributaries to a main glacier
global_tasks.compute_ref_t_stars Detects the best t* for the reference glaciers and writes them to disk
global_tasks.get_ref_mb_glaciers Get the list of glaciers we have valid mass balance measurements for.
global_tasks.write_centerlines_to_shape Write the centerlines in a shapefile.
global_tasks.compile_run_output Merge the output of the model runs of several gdirs into one file.
global_tasks.compile_climate_input Merge the climate input files in the glacier directories into one file.
global_tasks.compile_task_log Gathers the log output for the selected task(s)
global_tasks.compile_task_time Gathers the time needed for the selected task(s) to run
global_tasks.compile_glacier_statistics Gather as much statistics as possible about a list of glaciers.
global_tasks.compile_fixed_geometry_mass_balance Compiles a table of specific mass-balance timeseries for all glaciers.
global_tasks.compile_climate_statistics Gather as much statistics as possible about a list of glaciers.

Command line interface (CLI)

These commands are available:

  • oggm_netrc_credentials
  • oggm_prepro
  • oggm_benchmark
cli.prepro_levels.run_prepro_levels Generate the preprocessed OGGM glacier directories for this OGGM version

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 Geometry (line and widths) and flow rooting properties, but no thickness
Flowline A Centerline with additional properties: input to the FlowlineModel
MassBalanceModel Interface and common logic for all mass balance models used in OGGM.
FlowlineModel Interface to OGGM’s flowline models
FileModel Duck FlowlineModel which actually reads data out of a nc file.

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_directories('RGI60-11.00897',
   ...:                                           from_prepro_level=1,
   ...:                                           prepro_border=10)
   ...: 

In [8]: gdir = gdirs[0]  # init_glacier_directories 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]: 
['glacier_grid.json',
 'outlines.tar.gz',
 'log.txt',
 'diagnostics.json',
 'dem_source.txt',
 'dem.tif',
 'intersects.tar.gz']

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');
_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 = workflow.init_glacier_directories('RGI60-11.00897')[0]

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

For more information about how to use GlacierDirectories, visit our tutorial on the topic.

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.

catchments_intersects.shp
The intersections between catchments (shapefile) in the local map projection (Transverse Mercator).
centerlines.pkl
A list of oggm.Centerline instances, sorted by flow order.
climate_historical.nc
The historical monthly climate timeseries stored in a netCDF file.
climate_historical_daily.nc
The historical daily climate timeseries stored in a netCDF file.(only temperature is really changing on daily basis,precipitation is just assumed constant for every day
climate_info.json
Some information (dictionary) about the mass balance parameters for this glacier.
climate_monthly.nc
Deprecated: old name for climate_historical.
dem.tif
A geotiff file containing the DEM (reprojected into the local grid).This DEM is not smoothed or gap files, 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 downstream line geometry as well as the bed shape computed from a parabolic fit.
elevation_band_flowline.csv
A table containing the Huss&Farinotti 2012 squeezed flowlines.
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.
glacier_mask.tif
A glacier mask geotiff file with the same extend and projection as the dem.tif. This geotiff has value 1 at glaciated grid points and value 0 at unglaciated points.
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.
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 Interface and common logic for all mass balance models used in OGGM.
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.

Model Flowlines

Flowlines are what represent a glacier in OGGM. During the preprocessing stage of a glacier simulation these lines evolve from simple topography following lines to more abstract objects within the model.

This chapter provides an overview on the different line types, how to access the model flowlines and their attributes.

Hierarchy

First a short summary of the evolution or the hierarchy of these different flowline stages:

geometrical centerlines
These lines are calculated following the algorithm of Kienholz et al., (2014). They determine where the subsequent flowlines are situated and how many of them are needed for one glacier. The main task for this calculation is compute_centerlines(). These lines are not stored within a specific OGGM class but are stored within the GlacierDirectory as shapely objects.
inversion flowlines
To use the geometrical lines within the model they are transformed to Centerline objects. This is done within initialize_flowlines() where the lines are projected onto a xy-grid with regular spaced line-points. This step also takes care of tributary junctions. These lines are then in a later step also used for the bed inversion, hence their name.
downstream line
This line extends the glacier’s centerline along the unglaciated glacier bed which is necessary for advancing glaciers. This line is calculated along the valley floor to the end of the domain boundary within compute_downstream_line().
model flowlines
The dynamical OGGM model finally uses lines of the class Flowline for all simulations. These flowlines are created by init_present_time_glacier(), which combines information from several preprocessing steps, the downstream line and the bed inversion. From a user perspective, especially if preprocessed directories are used, these model flowlines are the most important ones and further informations on the class interface and attributes are given below.

Access

The easiest and most common way to access the model flowlines of a glacier is with its GlacierDirectory. For this we initialize a minimal working example including the model flowlines. This is achieved by choosing preprocessing level 4.

In [22]: import os

In [23]: import oggm

In [24]: from oggm import cfg, workflow

In [25]: from oggm.utils import gettempdir

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

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

In [28]: gdirs = workflow.init_glacier_directories('RGI60-11.00897',
   ....:                                           from_prepro_level=5,
   ....:                                           prepro_border=80)
   ....: 

In [29]: gdir = gdirs[0]  # init_glacier_directories always returns a list

In [30]: fls = gdir.read_pickle('model_flowlines')

In [31]: fls
Out[31]: 
[<oggm.core.flowline.MixedBedFlowline at 0x7f5c49f48438>,
 <oggm.core.flowline.MixedBedFlowline at 0x7f5c49f43f98>,
 <oggm.core.flowline.MixedBedFlowline at 0x7f5c49f48a90>]

In [32]: [fl.order for fl in fls]
Out[32]: [0, 0, 1]

This glacier has three flowlines of type MixedBedFlowline provided as a list. And the flowlines are ordered by ascending Strahler numbers, where the last entry in the list is always the longest and very often most relevant flowline of that glacier.

Type of model flowlines

Flowline A Centerline with additional properties: input to the FlowlineModel
MixedBedFlowline A Flowline which can take a combination of different shapes (default)
ParabolicBedFlowline A parabolic shaped Flowline with one degree of freedom
RectangularBedFlowline Simple shaped Flowline, glacier width does not change with ice thickness
TrapezoidalBedFlowline A Flowline with trapezoidal shape and two degrees of freedom

Important flowline functions

centerlines.initialize_flowlines Computes more physical Inversion Flowlines from geometrical Centerlines
centerlines.compute_centerlines Compute the centerlines following Kienholz et al., (2014).
centerlines.compute_downstream_line Computes the Flowline along the unglaciated downstream topography
flowline.init_present_time_glacier Merges data from preprocessing tasks.