Step-by-Step guide to building preprocessed directories from scratch#
This guide is for anyone who wants to learn how to start from scratch with building preprocessed directories, the same way OGGM does for its users. You might even use your own data. This tutorial assumes you are already familiar with OGGM basics. We won’t go into every little detail here, but you’ll find links for more information if you’re interested.
We’ve structured the guide into five main sections, each dedicated to a different level of preprocessing. At the beginning of each section, we’ll outline the tasks to be performed, the data we’ll use, and provide links to related tutorials. Additionally, at the end of each section, we’ll share a corresponding prepro_base_url
. This URL allows you to start directly at that level with everything pre-setup, bypassing the need to complete the earlier steps. Plus, in the tutorial storing glacier directories for later use, we show you how to save your work. This way, you don’t have to redo everything from the beginning every time (many steps only need to be done once).
Tags: advanced, glacier-directory, workflow
Tip: There’s a lot to learn here. If you’re curious about a specific function and want to know more, just add a question mark (?) right after it, and you’ll see more details.
# Example: Getting help on the Python function 'print'
print?
Set-up#
First, let’s get everything ready to go. Here’s how we’ll do it:
Import Functions: We’ll start by importing the functions we need. Depending on which preprocessed level you’re working with, you might not need all of them.
Initialize OGGM: Next, we’ll set up OGGM and choose where we want to save our work (defining the working directory).
Choose a Glacier: Lastly, we’ll pick one glacier to focus on as our example.
Remember, these steps are important no matter which level you’re starting from!
from oggm import cfg, utils, workflow, tasks, DEFAULT_BASE_URL
import geopandas as gpd
import numpy as np
import os
# we always need to initialzie and define a working directory
cfg.initialize(logging_level='WARNING')
cfg.PATHS['working_dir'] = utils.gettempdir(dirname='OGGM-full_prepro_elevation_bands', reset=True)
2024-08-25 21:00:04: oggm.cfg: Reading default parameters from the OGGM `params.cfg` configuration file.
2024-08-25 21:00:04: oggm.cfg: Multiprocessing switched OFF according to the parameter file.
2024-08-25 21:00:04: oggm.cfg: Multiprocessing: using all available processors (N=4)
# Our example glacier
rgi_ids = ['RGI60-11.00897'] # Hintereisferner
rgi_region = '11' # this must fit to example glacier(s), if starting from level 0
# This section is only for future developments of the tutorial (e.g. updateing for new OGGM releases)
# Test if prepro_base_url valid for both flowline_type_to_use, see level 2.
# In total four complete executions of the notebook:
# (load_from_prepro_base_url=False/True and flowline_type_to_use = 'elevation_band'/'centerline')
load_from_prepro_base_url = False
Level 0#
Tasks:
Define the
rgi_id
for your glacier directorygdir
.Define the map projection of the glacier directory
Add an outline of the glacier.
Optionally add intersects to other outlines.
Data used:
Glacier outline
Optionally intersects
Related Tutorials:
working_with_rgi will show you how to read glacier outline files and prepare them for a run
use_your_own_inventory: Use custom glacier inventories with OGGM
# load all RGI outlines for our region and extract the example glaciers
rgidf = gpd.read_file(utils.get_rgi_region_file(rgi_region, version='62'))
rgidf = rgidf[np.isin(rgidf.RGIId, rgi_ids)]
# We also take care of intersects for this RGI version
cfg.set_intersects_db(utils.get_rgi_intersects_region_file(rgi_region, version='62'))
# set the used projection used for gdir, options 'tmerc' or 'utm'
cfg.PARAMS['map_proj'] = cfg.PARAMS['map_proj'] # default is 'tmerc'
gdirs = workflow.init_glacier_directories(rgidf, reset=True, force=True)
2024-08-25 21:00:11: oggm.workflow: Execute entity tasks [GlacierDirectory] on 1 glaciers
# Instruction for beginning with existing OGGM's preprocessed directories
if load_from_prepro_base_url:
# to start from level 0 you can do
prepro_base_url_L0 = 'https://cluster.klima.uni-bremen.de/~oggm/gdirs/oggm_v1.6/L1-L2_files/elev_bands/'
gdirs = workflow.init_glacier_directories(rgi_ids,
from_prepro_level=0,
prepro_base_url=prepro_base_url_L0,
prepro_border=80, # could be 10, 80, 160 or 240
reset=True,
force=True,
)
Level 1#
Tasks:
Define the border around the outline.
Define the local grid resolution, which will also set the resolution for the flowlines.
Add the digital elevation model DEM.
Set up a local grid for each
gdir
.
Data used:
DEM file
Related Tutorials:
dem_sources: Create local topography maps from different DEM sources with OGGM
rgitopo_rgi6: RGI-TOPO for RGI v6.0
# define the border, we keep the default here
cfg.PARAMS['border'] = cfg.PARAMS['border']
# set the method for determining the local grid resolution
cfg.PARAMS['grid_dx_method'] = cfg.PARAMS['grid_dx_method'] # The default method is 'square', which determines the grid spacing (dx) based on the glacier's outline area.
cfg.PARAMS['fixed_dx'] = cfg.PARAMS['fixed_dx'] # This allows setting a specific resolution in meters. It's applicable only when grid_dx_method is set to 'fixed'.
# set the DEM source to use
source = 'COPDEM90' # we use COPDEM here
# this task adds the DEM and defines the local grid
workflow.execute_entity_task(tasks.define_glacier_region,
gdirs,
source=source);
2024-08-25 21:00:11: oggm.workflow: Execute entity tasks [define_glacier_region] on 1 glaciers
# Instruction for beginning with existing OGGM's preprocessed directories
if load_from_prepro_base_url:
# to start from level 1 you can do
prepro_base_url_L1 = 'https://cluster.klima.uni-bremen.de/~oggm/gdirs/oggm_v1.6/L1-L2_files/elev_bands/'
gdirs = workflow.init_glacier_directories(rgi_ids,
from_prepro_level=1,
prepro_base_url=prepro_base_url_L1,
prepro_border=80, # could be 10, 80, 160 or 240
reset=True,
force=True,
)
Level 2#
Tasks:
Choose the type of flowline to use.
Create the flowlines surface structure, including surface height and width.
Create the downstream flowline, which starts from the glacier’s terminus and extends downstream.
Optionally you can bring in extra data from the OGGM-shop and bin it to the elevation band flowline.
Data used:
Outline
DEM
Optional: additional datasets
Related Tutorials:
elevation_bands_vs_centerlines: Differences between “elevation band” and “centerline” flowlines
numeric_solvers: Understand the difference between the ice dynamic solvers in OGGM
oggm_shop: OGGM-Shop and Glacier Directories in OGGM
ingest_gridded_data_on_flowlines: ingest gridded products such as ice velocity into OGGM (and compare them with model output)
flowline_type_to_use = 'elevation_band' # you can also select 'centerline' here
if flowline_type_to_use == 'elevation_band':
elevation_band_task_list = [
tasks.simple_glacier_masks,
tasks.elevation_band_flowline,
tasks.fixed_dx_elevation_band_flowline,
tasks.compute_downstream_line,
tasks.compute_downstream_bedshape,
]
for task in elevation_band_task_list:
workflow.execute_entity_task(task, gdirs);
elif flowline_type_to_use == 'centerline':
# for centerline we can use parabola downstream line
cfg.PARAMS['downstream_line_shape'] = 'parabola'
centerline_task_list = [
tasks.glacier_masks,
tasks.compute_centerlines,
tasks.initialize_flowlines,
tasks.catchment_area,
tasks.catchment_intersections,
tasks.catchment_width_geom,
tasks.catchment_width_correction,
tasks.compute_downstream_line,
tasks.compute_downstream_bedshape,
]
for task in centerline_task_list:
workflow.execute_entity_task(task, gdirs);
else:
raise ValueError(f"Unknown flowline type '{flowline_type_to_use}'! Select 'elevation_band' or 'centerline'!")
2024-08-25 21:00:11: oggm.workflow: Execute entity tasks [simple_glacier_masks] on 1 glaciers
2024-08-25 21:00:11: oggm.workflow: Execute entity tasks [elevation_band_flowline] on 1 glaciers
2024-08-25 21:00:12: oggm.workflow: Execute entity tasks [fixed_dx_elevation_band_flowline] on 1 glaciers
2024-08-25 21:00:12: oggm.workflow: Execute entity tasks [compute_downstream_line] on 1 glaciers
2024-08-25 21:00:12: oggm.workflow: Execute entity tasks [compute_downstream_bedshape] on 1 glaciers
# Instruction for beginning with existing OGGM's preprocessed directories
if load_from_prepro_base_url:
# to start from level 2 we need to distinguish between the flowline types
if flowline_type_to_use == 'elevation_band':
prepro_base_url_L2 = 'https://cluster.klima.uni-bremen.de/~oggm/gdirs/oggm_v1.6/L1-L2_files/2023.2/elev_bands_w_data/'
elif flowline_type_to_use == 'centerline':
prepro_base_url_L2 = 'https://cluster.klima.uni-bremen.de/~oggm/gdirs/oggm_v1.6/L1-L2_files/centerlines/'
else:
raise ValueError(f"Unknown flowline type '{flowline_type_to_use}'! Select 'elevation_band' or 'centerline'!")
gdirs = workflow.init_glacier_directories(rgi_ids,
from_prepro_level=2,
prepro_base_url=prepro_base_url_L2,
prepro_border=80, # could be 10, 80, 160 or 240
reset=True,
force=True,
)
Level 3#
Tasks:
Add baseline climate data to gdir.
Calibrate the mass balance model statically (without considering glacier dynamics) using geodetic observations. This involves the calibration of
melt_f
,prcp_fac
andtemp_bias
.Conduct an inversion for the glacier’s bed topography. Including the calibration of
glen_a
andfs
by matching to the total volume estimate.Create the dynamic flowline for dynamic simulation runs.
Data used:
Climate data (default source: GSWP3_W5E5)
Geodetic mass balance observations (default source: Hugonnet al. 2021)
Volume estimate (default source: Farinotti et al. (2019))
Related Tutorials:
oggm_shop: OGGM-Shop and Glacier Directories in OGGM
massbalance_calibration: A look into the new mass balance calibration in OGGM v1.6
massbalance_perturbation: Mass balance parameter perturbation experiments
inversion: run the OGGM ice thickness inversion model with various ice parameters
# define the climate data to use, we keep the default
cfg.PARAMS['baseline_climate'] = cfg.PARAMS['baseline_climate']
# add climate data to gdir
workflow.execute_entity_task(tasks.process_climate_data, gdirs);
# the default mb calibration
workflow.execute_entity_task(tasks.mb_calibration_from_geodetic_mb,
gdirs,
informed_threestep=True, # only available for 'GSWP3_W5E5'
);
# glacier bed inversion
workflow.execute_entity_task(tasks.apparent_mb_from_any_mb, gdirs);
workflow.calibrate_inversion_from_consensus(
gdirs,
apply_fs_on_mismatch=True,
error_on_mismatch=True, # if you running many glaciers some might not work
filter_inversion_output=True, # this partly filters the overdeepening due to
# the equilibrium assumption for retreating glaciers (see. Figure 5 of Maussion et al. 2019)
volume_m3_reference=None, # here you could provide your own total volume estimate in m3
);
# finally create the dynamic flowlines
workflow.execute_entity_task(tasks.init_present_time_glacier, gdirs);
2024-08-25 21:00:12: oggm.workflow: Execute entity tasks [process_climate_data] on 1 glaciers
2024-08-25 21:00:12: oggm.workflow: Execute entity tasks [mb_calibration_from_geodetic_mb] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [apparent_mb_from_any_mb] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Applying global task calibrate_inversion_from_consensus on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Consensus estimate optimisation with A factor: 0.1 and fs: 0
2024-08-25 21:00:13: oggm.workflow: Applying global task inversion_tasks on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [prepare_for_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [mass_conservation_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [filter_inversion_output] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [get_inversion_volume] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Consensus estimate optimisation with A factor: 10.0 and fs: 0
2024-08-25 21:00:13: oggm.workflow: Applying global task inversion_tasks on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [prepare_for_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [mass_conservation_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [filter_inversion_output] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [get_inversion_volume] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Consensus estimate optimisation with A factor: 7.876802795318141 and fs: 0
2024-08-25 21:00:13: oggm.workflow: Applying global task inversion_tasks on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [prepare_for_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [mass_conservation_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [filter_inversion_output] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [get_inversion_volume] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Consensus estimate optimisation with A factor: 3.98840139765907 and fs: 0
2024-08-25 21:00:13: oggm.workflow: Applying global task inversion_tasks on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [prepare_for_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [mass_conservation_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [filter_inversion_output] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [get_inversion_volume] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Consensus estimate optimisation with A factor: 2.044200698829535 and fs: 0
2024-08-25 21:00:13: oggm.workflow: Applying global task inversion_tasks on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [prepare_for_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [mass_conservation_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [filter_inversion_output] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [get_inversion_volume] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Consensus estimate optimisation with A factor: 2.869309409059597 and fs: 0
2024-08-25 21:00:13: oggm.workflow: Applying global task inversion_tasks on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [prepare_for_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [mass_conservation_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [filter_inversion_output] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [get_inversion_volume] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Consensus estimate optimisation with A factor: 2.7157333579518252 and fs: 0
2024-08-25 21:00:13: oggm.workflow: Applying global task inversion_tasks on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [prepare_for_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [mass_conservation_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [filter_inversion_output] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [get_inversion_volume] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Consensus estimate optimisation with A factor: 2.6897121901531937 and fs: 0
2024-08-25 21:00:13: oggm.workflow: Applying global task inversion_tasks on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [prepare_for_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [mass_conservation_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [filter_inversion_output] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [get_inversion_volume] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: calibrate_inversion_from_consensus converged after 7 iterations and fs=0. The resulting Glen A factor is 2.6897121901531937.
2024-08-25 21:00:13: oggm.workflow: Applying global task inversion_tasks on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [prepare_for_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [mass_conservation_inversion] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [filter_inversion_output] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [get_inversion_volume] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [init_present_time_glacier] on 1 glaciers
Guidance on utilizing various baseline climates:
Currently, OGGM supports a variety of baseline climates, including ‘CRU’, ‘HISTALP’, ‘W5E5’, ‘GSWP3_W5E5’ (the default), ‘ERA5’, ‘ERA5L’, ‘CERA’, ‘ERA5dr’, and ‘ERA5L-HMA’. Although switching between these datasets is straightforward, calibrating the mass balance model according to each dataset is more complex. For instance, you’ll need to choose a default precipitation factor that suits both your selected climate dataset and your specific region. Additionally, you must determine the best method to calibrate the mass balance parameters. For a comprehensive guide on the available options, explanations, and how to incorporate your own geodetic observations, please refer to the tutorial massbalance_calibration.
Here’s an example of using the ERA5 dataset:
# define the baseline climate and add it
cfg.PARAMS['baseline_climate'] = 'ERA5'
workflow.execute_entity_task(tasks.process_climate_data, gdirs);
# define the default precipitation factor
cfg.PARAMS['prcp_fac'] = 1.6 # Note: This is not a universial value!
cfg.PARAMS['use_winter_prcp_fac'] = False # This option is only available for 'GSWP3_W5E5'
cfg.PARAMS['use_temp_bias_from_file'] = False # This option is only available for 'GSWP3_W5E5'
# an example of static calibration for mass balance, more options are available in the tutorial
workflow.execute_entity_task(tasks.mb_calibration_from_geodetic_mb,
gdirs,
calibrate_param1='melt_f',
calibrate_param2='prcp_fac',
calibrate_param3='temp_bias')
You can also utilize your own climate data. However, you will need to either convert your data into a specific format (for an example, see OGGM/oggm-sample-data ->test-files/histalp_merged_hef.nc) or create your own tasks.process_climate_data
function. Here’s how you might do this:
cfg.PARAMS['baseline_climate'] = 'CUSTOM'
cfg.PATHS['climate_file'] = path_to_the_climate_file
workflow.execute_entity_task(tasks.process_climate_data, gdirs);
# proceed with defining the default precipitation factor and mass balance calibration as shown above
# Instruction for beginning with existing OGGM's preprocessed directories
if load_from_prepro_base_url:
# to start from level 3 you can do
if flowline_type_to_use == 'elevation_band':
prepro_base_url_L3 = 'https://cluster.klima.uni-bremen.de/~oggm/gdirs/oggm_v1.6/L3-L5_files/2023.3/elev_bands/W5E5/'
elif flowline_type_to_use == 'centerline':
prepro_base_url_L3 = 'https://cluster.klima.uni-bremen.de/~oggm/gdirs/oggm_v1.6/L3-L5_files/2023.3/centerlines/W5E5/'
else:
raise ValueError(f"Unknown flowline type '{flowline_type_to_use}'! Select 'elevation_band' or 'centerline'!")
gdirs = workflow.init_glacier_directories(rgi_ids,
from_prepro_level=3,
prepro_base_url=prepro_base_url_L3,
prepro_border=80, # could be 80 or 160
reset=True,
force=True,
)
Level 4#
Tasks:
Initialize the current state of the glacier without a dynamic spinup. This method, default until version 1.6., is mainly for comparison purposes and can often be skipped.
Initialize the current glacier state with a dynamic spinup. This process includes a dynamic calibration of the mass balance. It’s important to note that this option isn’t available for centerlines in the current OGGM preprocessed directories, meaning it hasn’t been tested or analyzed.
Data used:
Geodetic mass balance observations (default source: Hugonnet al. 2021)
Related Tutorials:
10 minutes to… the new dynamical spinup in OGGM v1.6
dynamical_spinup: a deeper dive into the dynamical spinup for past simulations
numeric_solvers: Understand the difference between the ice dynamic solvers in OGGM
# set the ice dynamic solver depending on the flowline-type
if flowline_type_to_use == 'elevation_band':
cfg.PARAMS['evolution_model'] = 'SemiImplicit'
elif flowline_type_to_use == 'centerline':
cfg.PARAMS['evolution_model'] = 'FluxBased'
else:
raise ValueError(f"Unknown flowline type '{flowline_type_to_use}'! Select 'elevation_band' or 'centerline'!")
# get the start and end year of the selected baseline
y0 = gdirs[0].get_climate_info()['baseline_yr_0']
ye = gdirs[0].get_climate_info()['baseline_yr_1'] + 1 # run really to the end until 1.1.
# 'static' initialisation
workflow.execute_entity_task(tasks.run_from_climate_data, gdirs,
min_ys=y0, ye=ye,
fixed_geometry_spinup_yr=None, # here you could add a static spinup if you want
output_filesuffix='_historical')
# 'dynamic' initialisation, including dynamic mb calibration
dynamic_spinup_start_year = 1979
minimise_for = 'area' # other option would be 'volume'
workflow.execute_entity_task(
tasks.run_dynamic_melt_f_calibration, gdirs,
err_dmdtda_scaling_factor=0.2, # by default we reduce the mass balance error for accounting for
# corrleated uncertainties on a regional scale
ys=dynamic_spinup_start_year, ye=ye,
kwargs_run_function={'minimise_for': minimise_for},
ignore_errors=True,
kwargs_fallback_function={'minimise_for': minimise_for},
output_filesuffix='_spinup_historical',
);
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [run_from_climate_data] on 1 glaciers
2024-08-25 21:00:13: oggm.workflow: Execute entity tasks [run_dynamic_melt_f_calibration] on 1 glaciers
# Instruction for beginning with existing OGGM's preprocessed directories
if load_from_prepro_base_url:
# to start from level 4 you can do
if flowline_type_to_use == 'elevation_band':
prepro_base_url_L4 = DEFAULT_BASE_URL
elif flowline_type_to_use == 'centerline':
prepro_base_url_L4 = 'https://cluster.klima.uni-bremen.de/~oggm/gdirs/oggm_v1.6/L3-L5_files/2023.3/centerlines/W5E5/'
else:
raise ValueError(f"Unknown flowline type '{flowline_type_to_use}'! Select 'elevation_band' or 'centerline'!")
gdirs = workflow.init_glacier_directories(rgi_ids,
from_prepro_level=4,
prepro_base_url=prepro_base_url_L4,
prepro_border=80, # could be 80 or 160
reset=True,
force=True,
)
Level 5#
Tasks:
Retain only the data necessary for future projection runs to conserve disk space. At this stage, it’s not possible to revisit the preprocessing steps from earlier levels, but all required information for conducting future projection runs is preserved.
Data used:
No additional data is needed for this level.
Related Tutorials:
store_and_compress_glacierdirs: storing glacier directories for later use
mini_base_dir = os.path.join(cfg.PATHS['working_dir'],
'mini_per_glacier')
mini_gdirs = workflow.execute_entity_task(tasks.copy_to_basedir, gdirs,
base_dir=mini_base_dir,
setup='run/spinup')
2024-08-25 21:00:16: oggm.workflow: Execute entity tasks [copy_to_basedir] on 1 glaciers
# Instruction for beginning with existing OGGM's preprocessed directories
if load_from_prepro_base_url:
# to start from level 5 you can do
if flowline_type_to_use == 'elevation_band':
prepro_base_url_L5 = DEFAULT_BASE_URL
elif flowline_type_to_use == 'centerline':
prepro_base_url_L5 = 'https://cluster.klima.uni-bremen.de/~oggm/gdirs/oggm_v1.6/L3-L5_files/2023.3/centerlines/W5E5/'
else:
raise ValueError(f"Unknown flowline type '{flowline_type_to_use}'! Select 'elevation_band' or 'centerline'!")
gdirs = workflow.init_glacier_directories(rgi_ids,
from_prepro_level=5,
prepro_base_url=prepro_base_url_L5,
prepro_border=80, # could be 80 or 160
reset=True,
force=True,
)
And that’s it! We’ve successfully recreated all the preprocessed levels offered by OGGM. Remember, if you prefer, you can bypass all previous steps and jump straight into your future projections from Level 5. Happy modeling!
What’s next?#
return to the OGGM documentation
back to the table of contents