pygem package#

Subpackages#

Submodules#

pygem.class_climate module#

Python Glacier Evolution Model (PyGEM)

copyright © 2018 David Rounce <drounce@cmu.edu

Distrubted under the MIT lisence

class of climate data and functions associated with manipulating the dataset to be in the proper format

class pygem.class_climate.GCM(name='', sim_climate_scenario='', realization=None)[source]#

Bases: object

Global climate model data properties and functions used to automatically retrieve data.

Attributes:
namestr

name of climate dataset.

sim_climate_scenariostr

rcp or ssp sim_climate_scenario (example: ‘rcp26’ or ‘ssp585’)

realizationstr

realization from large ensemble (example: ‘1011.001’ or ‘1301.020’)

Methods

importGCMfxnearestneighbor_xarray(filename, ...)

Import time invariant (constant) variables and extract nearest neighbor.

importGCMvarnearestneighbor_xarray(filename, ...)

Import time series of variables and extract nearest neighbor.

importGCMfxnearestneighbor_xarray(filename, vn, main_glac_rgi)[source]#

Import time invariant (constant) variables and extract nearest neighbor.

Note: cmip5 data used surface height, while ERA-Interim data is geopotential

Parameters:
filenamestr

filename of variable

variablenamestr

variable name

main_glac_rgipandas dataframe

dataframe containing relevant rgi glacier information

Returns:
glac_variablenumpy array

array of nearest neighbor values for all the glaciers in model run (rows=glaciers, column=variable)

importGCMvarnearestneighbor_xarray(filename, vn, main_glac_rgi, dates_table, realizations=['r1i1p1f1', 'r4i1p1f1'])[source]#

Import time series of variables and extract nearest neighbor.

Note: “NG” refers to a homogenized “new generation” of products from ETH-Zurich.

The function is setup to select netcdf data using the dimensions: time, latitude, longitude (in that order). Prior to running the script, the user must check that this is the correct order of the dimensions and the user should open the netcdf file to determine the names of each dimension as they may vary.

Parameters:
filenamestr

filename of variable

vnstr

variable name

main_glac_rgipandas dataframe

dataframe containing relevant rgi glacier information

dates_table: pandas dataframe

dataframe containing dates of model run

Returns:
glac_variable_seriesnumpy array

array of nearest neighbor values for all the glaciers in model run (rows=glaciers, columns=time series)

time_seriesnumpy array

array of dates associated with the meteorological data (may differ slightly from those in table for monthly timestep, i.e., be from the beginning/middle/end of month)

pygem.gcmbiasadj module#

Python Glacier Evolution Model (PyGEM)

copyright © 2018 David Rounce <drounce@cmu.edu>

Distrubted under the MIT lisence

Run bias adjustments a given climate dataset

pygem.gcmbiasadj.annual_avg_2darray(x)[source]#

Annual average of dataset, where columns are a monthly timeseries (temperature)

pygem.gcmbiasadj.annual_sum_2darray(x)[source]#

Annual sum of dataset, where columns are a monthly timeseries (precipitation)

pygem.gcmbiasadj.monthly_avg_2darray(x)[source]#

Monthly average for a given 2d dataset where columns are monthly timeseries

pygem.gcmbiasadj.monthly_avg_array_rolled(ref_array, dates_table_ref, dates_table, sim_startyear, ref_startyear)[source]#

Monthly average array from reference data rolled to ensure proper months

Parameters:
ref_arraynp.array

time series of reference lapse rates

dates_table_refpd.DataFrame

dates table for reference time period

dates_tablepd.DataFrame

dates_table for GCM time period

Returns:
gcm_arraynp.array

gcm climate data based on monthly average of reference data

pygem.gcmbiasadj.monthly_std_2darray(x)[source]#

Monthly standard deviation for a given 2d dataset where columns are monthly timeseries

pygem.gcmbiasadj.prec_biasadj_HH2015(ref_prec, ref_elev, gcm_prec, dates_table_ref, dates_table, sim_startyear, ref_startyear)[source]#

Huss and Hock (2015) precipitation bias correction based on mean (multiplicative)

Parameters:
ref_precnp.array

time series of reference precipitation

gcm_precnp.array

time series of GCM precipitation

dates_table_refpd.DataFrame

dates table for reference time period

dates_tablepd.DataFrame

dates_table for GCM time period

gcm_elev_biasadjfloat

new gcm elevation is the elevation of the reference climate dataset

Returns:
gcm_prec_biasadjnp.array

GCM precipitation bias corrected to the reference climate dataset according to Huss and Hock (2015)

pygem.gcmbiasadj.prec_biasadj_QDM(ref_prec, ref_elev, gcm_prec, dates_table_ref, dates_table, sim_startyear, ref_startyear)[source]#
Cannon et al. (2015) precipitation bias correction based on quantile delta mapping

Also see Lader et al. (2017) another use case

Perform a quantile delta mapping bias-correction procedure on precipitation.

This function operates by multiplying reference precipitation by a ratio of

the projected and future gcm precipitations at the same percentiles (e.g., ref_prec * gcm_projected/gcm_historic, with all values at same percentile).

Quantile delta mapping is generally viewed as more capable of capturing

climatic extemes at the lowest and highest quantiles (e.g., 0.01% and 99.9%) compared to standard quantile mapping (which constructs a transfer function using only reference and historic climate data, requiring extrapolations for projected values lying outside the reference and historic datasets). See Cannon et al. (2015) Sections 2 and 3 for further explanation.

Parameters:
ref_precpandas dataframe

dataframe containing reference climate precipitation data

gcm_precnp.array

time series of GCM precipitation

dates_table_refpd.DataFrame

dates table for reference time period

dates_tablepd.DataFrame

dates_table for GCM time period

Returns:
gcm_prec_biasadjnumpy ndarray

ndarray that contains bias-corrected future gcm precipitation data

gcm_elev_biasadjfloat

new gcm elevation is the elevation of the reference climate dataset

pygem.gcmbiasadj.prec_biasadj_opt1(ref_prec, ref_elev, gcm_prec, dates_table_ref, dates_table, sim_startyear, ref_startyear)[source]#

Precipitation bias correction based on mean with limited maximum

Parameters:
ref_precnp.array

time series of reference precipitation

gcm_precnp.array

time series of GCM precipitation

dates_table_refpd.DataFrame

dates table for reference time period

dates_tablepd.DataFrame

dates_table for GCM time period

Returns:
gcm_prec_biasadjnp.array

GCM precipitation bias corrected to the reference climate dataset according to Huss and Hock (2015)

gcm_elev_biasadjfloat

new gcm elevation is the elevation of the reference climate dataset

pygem.gcmbiasadj.temp_biasadj_HH2015(ref_temp, ref_elev, gcm_temp, dates_table_ref, dates_table, sim_startyear, ref_startyear, debug=False)[source]#

Huss and Hock (2015) temperature bias correction based on mean and interannual variability

Note: the mean over the reference period will only equal the mean of the gcm for the same time period when the GCM time series is run for the same period, i.e., due to the 25-year moving average, the mean gcm temps from 2000-2019 will differ if using a reference period of 2000-2020 to bias adjust gcm temps from 2000-2100.

Parameters:
ref_tempnp.array

time series of reference temperature

gcm_tempnp.array

time series of GCM temperature

dates_table_refpd.DataFrame

dates table for reference time period

dates_tablepd.DataFrame

dates_table for GCM time period

Returns:
gcm_temp_biasadjnp.array

GCM temperature bias corrected to the reference climate dataset according to Huss and Hock (2015)

gcm_elev_biasadjfloat

new gcm elevation is the elevation of the reference climate dataset

pygem.gcmbiasadj.temp_biasadj_QDM(ref_temp, ref_elev, gcm_temp, dates_table_ref, dates_table, sim_startyear, ref_startyear)[source]#
Cannon et al. (2015) temperature bias correction based on quantile delta mapping

Also see Lader et al. (2017) for further documentation

Perform a quantile delta mapping bias-correction procedure on temperature.

This function operates by multiplying reference temperature by a ratio of

the projected and future gcm temperature at the same percentiles (e.g., ref_temp * gcm_projected/gcm_historic, with all values at same percentile).

Quantile delta mapping is generally viewed as more capable of capturing

climatic extemes at the lowest and highest quantiles (e.g., 0.01% and 99.9%) compared to standard quantile mapping (which constructs a transfer function using only reference and historic climate data, requiring extrapolations for projected values lying outside the reference and historic datasets). See Cannon et al. (2015) Sections 2 and 3 for further explanation.

Parameters:
ref_temppandas dataframe

dataframe containing reference climate temperature data

gcm_tempnp.array

time series of GCM temperature

dates_table_refpd.DataFrame

dates table for reference time period

dates_tablepd.DataFrame

dates_table for GCM time period

Returns:
gcm_temp_biasadjnumpy ndarray

ndarray that contains bias-corrected future gcm temperature data

gcm_elev_biasadjfloat

new gcm elevation is the elevation of the reference climate dataset

pygem.glacierdynamics module#

Python Glacier Evolution Model (PyGEM)

copyright © 2018 David Rounce <drounce@cmu.edu

Distrubted under the MIT lisence

class pygem.glacierdynamics.MassRedistributionCurveModel(flowlines, mb_model=None, y0=0.0, glen_a=None, fs=0.0, is_tidewater=False, water_level=0, inplace=False, debug=True, option_areaconstant=False, spinupyears=0, constantarea_years=0, **kwargs)[source]#

Bases: FlowlineModel

Glacier geometry updated using mass redistribution curves; also known as the “delta-h method”

This uses mass redistribution curves from Huss et al. (2010) to update the glacier geometry

Attributes:
area_km2
area_m2
length_m
mb_model
volume_bsl_km3
volume_bsl_m3
volume_bwl_km3
volume_bwl_m3
volume_km3
volume_m3
yr

Methods

check_domain_end()

Returns False if the glacier reaches the domains bound.

get_mb(heights[, year, fl_id, fls])

Get the mass balance at the requested height and time.

reset_flowlines(flowlines[, inplace, ...])

Reset the initial model flowlines

reset_y0(y0)

Reset the initial model time

run_until(y1[, run_single_year])

Runs the model from the current year up to a given year date y1.

run_until_and_store(y1[, run_path, ...])

Runs the model and returns intermediate steps in xarray datasets.

run_until_equilibrium([rate, ystep, max_ite])

Runs the model until an equilibrium state is reached.

step(dt)

Advance the numerical simulation of one single step.

to_geometry_netcdf(path)

Creates a netcdf group file storing the state of the model.

updategeometry(year[, debug])

Update geometry for a given year

run_until(y1, run_single_year=False)[source]#

Runs the model from the current year up to a given year date y1.

This function runs the model for the time difference y1-self.y0 If self.y0 has not been specified at some point, it is 0 and y1 will be the time span in years to run the model for.

Parameters:
y1float

Upper time span for how long the model should run

run_until_and_store(y1, run_path=None, diag_path=None, store_monthly_step=None)[source]#

Runs the model and returns intermediate steps in xarray datasets.

This function repeatedly calls FlowlineModel.run_until for either monthly or yearly time steps up till the upper time boundary y1.

Parameters:
y1int

Upper time span for how long the model should run (needs to be a full year)

run_pathstr

Path and filename where to store the model run dataset

diag_pathstr

Path and filename where to store the model diagnostics dataset

store_monthly_stepBool

If True (False) model diagnostics will be stored monthly (yearly). If unspecified, we follow the update of the MB model, which defaults to yearly (see __init__).

Returns:
run_dsxarray.Dataset

stores the entire glacier geometry. It is useful to visualize the glacier geometry or to restart a new run from a modelled geometry. The glacier state is stored at the begining of each hydrological year (not in between in order to spare disk space).

diag_dsxarray.Dataset

stores a few diagnostic variables such as the volume, area, length and ELA of the glacier.

updategeometry(year, debug=False)[source]#

Update geometry for a given year

pygem.massbalance module#

Python Glacier Evolution Model (PyGEM)

copyright © 2018 David Rounce <drounce@cmu.edu

Distrubted under the MIT lisence

class pygem.massbalance.PyGEMMassBalance(gdir, modelprms, glacier_rgi_table, option_areaconstant=False, frontalablation_k=None, debug=False, debug_refreeze=False, fls=None, fl_id=0, heights=None, inversion_filter=False, ignore_debris=False)[source]#

Bases: MassBalanceModel

Mass-balance computed from the Python Glacier Evolution Model.

This mass balance accounts for ablation, accumulation, and refreezing.

This class implements the MassBalanceModel interface so that the dynamical model can use it.

Methods

ensure_mass_conservation(diag)

Ensure mass conservation that may result from using OGGM's glacier dynamics model.

get_annual_mb(heights[, year, fls, fl_id, ...])

FIXED FORMAT FOR THE FLOWLINE MODEL

get_annual_specific_mass_balance(fls, year)

Get annual specific mass balance from multiple flowlines.

get_ela([year])

Get the equilibrium line altitude for a given year.

get_monthly_mb(heights[, year, fl_id, fls])

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

get_specific_mb([heights, widths, fls, year])

Specific mass balance for a given glacier geometry.

get_year_index(year)

Get the index of a specified year in the list of years available.

is_year_valid(year)

Checks if a given date year be simulated by this model.

reset_state()

Reset any internal state of the model.

ensure_mass_conservation(diag)[source]#

Ensure mass conservation that may result from using OGGM’s glacier dynamics model. This will be resolved on an annual basis, and since the glacier dynamics are updated annually, the melt and runoff will be adjusted on a monthly-scale based on percent changes.

OGGM’s dynamic model limits mass loss based on the ice thickness and flux divergence. As a result, the actual volume change, glacier runoff, glacier melt, etc. may be less than that recorded by the mb_model. For PyGEM this is important because the glacier runoff and all parameters should be mass conserving.

Note: other dynamical models (e.g., mass redistribution curves, volume-length-area scaling) are based on the total volume change and therefore do not impose limitations like this because they do not estimate the flux divergence. As a result, they may systematically overestimate mass loss compared to OGGM’s dynamical model.

get_annual_mb(heights, year=None, fls=None, fl_id=None, debug=False, option_areaconstant=False)[source]#

FIXED FORMAT FOR THE FLOWLINE MODEL

Returns annual climatic mass balance [m ice per second]

Parameters:
heightsnp.array

elevation bins

yearint

year to get the annual mass balance for

Returns:
mbnp.array

mass balance for each bin [m ice per second]

get_year_index(year)[source]#

Get the index of a specified year in the list of years available.

class pygem.massbalance.PyGEMMassBalance_wrapper(gdir, mb_model_class=<class 'pygem.massbalance.PyGEMMassBalance'>, fl_str='inversion_flowlines', filesuffix='', **kwargs)[source]#

Bases: MassBalanceModel

Attributes:
mbmod

Methods

get_annual_mb(heights[, year, fls, fl_id, ...])

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

get_annual_specific_mass_balance(fls, year)

Get annual specific mass balance from multiple flowlines.

get_ela([year])

Get the equilibrium line altitude for a given year.

get_monthly_mb(heights[, year, fl_id, fls])

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

get_specific_mb([heights, widths, fls, year])

Specific mass balance for a given glacier geometry.

is_year_valid(year)

Checks if a given date year be simulated by this model.

reset_state()

Reset any internal state of the model.

get_annual_mb(heights, year=None, fls=None, fl_id=None, debug=True, option_areaconstant=False)[source]#

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

For some simpler mass balance models get_monthly_mb()` and `get_annual_mb() can be equivalent.

Units: [m s-1], or meters of ice per second

Note: year is optional because some simpler models have no time component.

Parameters:
heights: ndarray

the altitudes at which the mass balance will be computed

year: float, optional

the time (in the “floating year” convention)

fl_id: float, optional

the index of the flowline in the fls array (might be ignored by some MB models)

fls: list of flowline instances, optional

the flowlines array, in case the MB model implementation needs to know details about the glacier geometry at the moment the MB model is called

Returns:
the mass balance (same dim as heights) (units: [m s-1])
property mbmod#

pygem.oggm_compat module#

Python Glacier Evolution Model (PyGEM)

copyright © 2018 David Rounce <drounce@cmu.edu

Distrubted under the MIT lisence

PYGEM-OGGGM COMPATIBILITY FUNCTIONS

class pygem.oggm_compat.CompatGlacDir(rgiid)[source]#

Bases: object

class pygem.oggm_compat.RandomLinearMassBalance(gdir, grad=3.0, h_perc=60, sigma_ela=100.0, seed=None)[source]#

Bases: MassBalanceModel

Mass-balance as a linear function of altitude with random ELA.

This is a dummy MB model to illustrate how to program one.

The reference ELA is taken at a percentile altitude of the glacier. It then varies randomly from year to year.

This class implements the MassBalanceModel interface so that the dynamical model can use it. Even if you are not familiar with object oriented programming, I hope that the example below is simple enough.

Methods

get_annual_mb(heights[, year, fl_id])

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

get_annual_specific_mass_balance(fls, year)

Get annual specific mass balance from multiple flowlines.

get_ela([year])

Get the equilibrium line altitude for a given year.

get_monthly_mb(heights[, year, fl_id, fls])

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

get_random_ela_h(year)

This generates a random ELA for the requested year.

get_specific_mb([heights, widths, fls, year])

Specific mass balance for a given glacier geometry.

is_year_valid(year)

Checks if a given date year be simulated by this model.

reset_state()

Reset any internal state of the model.

get_annual_mb(heights, year=None, fl_id=None)[source]#

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

For some simpler mass balance models get_monthly_mb()` and `get_annual_mb() can be equivalent.

Units: [m s-1], or meters of ice per second

Note: year is optional because some simpler models have no time component.

Parameters:
heights: ndarray

the altitudes at which the mass balance will be computed

year: float, optional

the time (in the “floating year” convention)

fl_id: float, optional

the index of the flowline in the fls array (might be ignored by some MB models)

fls: list of flowline instances, optional

the flowlines array, in case the MB model implementation needs to know details about the glacier geometry at the moment the MB model is called

Returns:
the mass balance (same dim as heights) (units: [m s-1])
get_random_ela_h(year)[source]#

This generates a random ELA for the requested year.

Since we do not know which years are going to be asked for we generate them on the go.

pygem.oggm_compat.create_empty_glacier_directory(rgi_id)[source]#

Create empty GlacierDirectory for PyGEM’s alternative ice thickness products

Parameters:
rgi_idstr

the rgi id of the glacier (RGIv60-)

Returns:
a GlacierDirectory object
pygem.oggm_compat.get_glacier_zwh(gdir)[source]#

Computes this glaciers altitude, width and ice thickness.

Parameters:
gdirGlacierDirectory

the glacier to compute

Returns:
a dataframe with the requested data
pygem.oggm_compat.single_flowline_glacier_directory(rgi_id, reset=False, prepro_border=80, logging_level='WORKFLOW', has_internet=True, working_dir='/path/to/pygem/data///oggm_gdirs/')[source]#

Prepare a GlacierDirectory for PyGEM (single flowline to start with)

Parameters:
rgi_idstr

the rgi id of the glacier (RGIv60-)

resetbool

set to true to delete any pre-existing files. If false (the default), the directory won’t be re-downloaded if already available locally in order to spare time.

prepro_borderint

the size of the glacier map: 10, 80, 160, 240

Returns:
a GlacierDirectory object
pygem.oggm_compat.single_flowline_glacier_directory_with_calving(rgi_id, reset=False, prepro_border=80, k_calving=1, logging_level='WORKFLOW', has_internet=True, working_dir='/path/to/pygem/data//oggm_gdirs/', facorrected=True)[source]#

Prepare a GlacierDirectory for PyGEM (single flowline to start with)

k_calving is free variable!

Parameters:
rgi_idstr

the rgi id of the glacier

resetbool

set to true to delete any pre-existing files. If false (the default), the directory won’t be re-downloaded if already available locally in order to spare time.

prepro_borderint

the size of the glacier map: 10, 80, 160, 250

Returns
——-
a GlacierDirectory object
pygem.oggm_compat.update_cfg(updates, dict_name='PARAMS')[source]#

Update keys in the OGGMs config.

Parameters: dict (str): The dictionary in the config to update. updates (dict): Key-Value pairs to be updated.

Returns: None: The function updates cfg in place.

pygem.pygem_modelsetup module#

Python Glacier Evolution Model (PyGEM)

copyright © 2018 David Rounce <drounce@cmu.edu

Distrubted under the MIT lisence

List of functions used to set up different aspects of the model

pygem.pygem_modelsetup.datesmodelrun(startyear=2000, endyear=2019, spinupyears=0, option_wateryear='calendar')[source]#

Create table of year, month, day, water year, season and number of days in the month.

Parameters:
startyearint

starting year

endyearint

ending year

spinupyearsint

number of spinup years

Returns:
dates_tablepd.DataFrame

table where each row is a timestep and each column is attributes (date, year, wateryear, etc.) of that timestep

pygem.pygem_modelsetup.daysinmonth(year, month)[source]#

Return days in month based on the month and year

Parameters:
yearstr
monthstr
Returns:
integer of the days in the month
pygem.pygem_modelsetup.hypsometrystats(hyps_table, thickness_table)[source]#

Calculate the volume and mean associated with the hypsometry data.

Output is a series of the glacier volume [km**3] and mean elevation values [m a.s.l.].

pygem.pygem_modelsetup.import_Husstable(rgi_table, filepath, filedict, drop_col_names, indexname='GlacNo', option_shift_elevbins_20m=True)[source]#

Use the dictionary specified by the user to extract the desired variable. The files must be in the proper units (ice thickness [m], area [km2], width [km]) and should be pre-processed.

Output is a Pandas DataFrame of the variable for all the glaciers in the model run (rows = GlacNo, columns = elevation bins).

Line Profiling: Loading in the table takes the most time (~2.3 s)

pygem.pygem_modelsetup.selectglaciersrgitable(glac_no=None, rgi_regionsO1=None, rgi_regionsO2='all', rgi_glac_number='all', rgi_fp='/path/to/pygem/data//RGI/rgi60/00_rgi60_attribs/', rgi_cols_drop=['GLIMSId', 'BgnDate', 'EndDate', 'Status', 'Linkages', 'Name'], rgi_O1Id_colname='glacno', rgi_glacno_float_colname='RGIId_float', indexname='GlacNo', include_landterm=True, include_laketerm=True, include_tidewater=True, glac_no_skip=None, min_glac_area_km2=0, debug=False)[source]#

Select all glaciers to be used in the model run according to the regions and glacier numbers defined by the RGI glacier inventory. This function returns the rgi table associated with all of these glaciers.

glac_nolist of strings

list of strings of RGI glacier numbers (e.g., [‘1.00001’, ‘13.00001’])

rgi_regionsO1list of integers

list of integers of RGI order 1 regions (e.g., [1, 13])

rgi_regionsO2list of integers or ‘all’

list of integers of RGI order 2 regions or simply ‘all’ for all the order 2 regions

rgi_glac_numberlist of strings

list of RGI glacier numbers without the region (e.g., [‘00001’, ‘00002’])

Output: Pandas DataFrame of the glacier statistics for each glacier in the model run (rows = GlacNo, columns = glacier statistics)

pygem.pygem_modelsetup.split_list(lst, n=1, option_ordered=1, group_thousands=False)[source]#

Split list into batches for the supercomputer.

Parameters:
lstlist

List that you want to split into separate batches

nint

Number of batches to split glaciers into.

Returns:
lst_batcheslist

list of n lists that have sequential values in each list

Module contents#

Python Glacier Evolution Model (PyGEM)

copyright © 2018 David Rounce <drounce@cmu.edu

Distrubted under the MIT lisence