MimiBRICK.jl

MimiBRICK.jl is an implementation of the Building Blocks for Relevant Ice and Climate Knowledge (BRICK) semi-empirical model for sea-level change in the Mimi integrated modeling framework (https://www.mimiframework.org/). The Mimi modeling framework is a coding platform that facilitates coupling models and running coupled modeling experiments. MimiBRICK.jl is flexible, efficient, and modular, to facilitate incorporating BRICK into coupled models and integrated assessments of climate impacts in a modular fashion to provide global average as well as local sea-level projections. This focus on tight model coupling and integrated modeling is a key feature of MimiBRICK.jl and broader Mimi modeling framework.

This implementation includes examples for using observational data to calibrate the model, as well as various configurations in which MimiBRICK.jl is coupled to other climate model components. For users who do not wish to re-run computationally intensive model calibration algorithms, this implementation also includes scripts for using existing calibration output for standard future climate change scenarios, and examples downscaling these global projections for assessments of local impacts. Pre-run model calibration and simulation output can be found in the accompanying Zenodo repository.

MimiBRICK.calculate_trendsMethod
calculate_trends(model_output::Array{Float64,1}, obs_trends::DataFrame, start_year::Int, end_year::Int)

Calculate modeled sea level trends.

Description: This function calculates the trend in different modeled sea level contributions over a user-specified time period.

Function Arguments:

  model_output = A vector of modeled annual sea level rise values.
  obs_trends   = The matrix of observed sea level trends provided by the "load_calibration_data" function defined above.
  start_year   = First year of the model run.
  end_year     = Last year of the model run.
source
MimiBRICK.construct_brick_log_posteriorMethod
construct_brick_log_posterior(f_run_model!; model_start_year::Int=1850, calibration_end_year::Int=2017, joint_antarctic_prior::Bool=false)

Calculate log posterior for brick.

Description: This creates a function that calculates the log-posterior probability of the uncertain model, initial condition, and statistical process parameters.

Function Arguments:

f_run_model           = A function that runs the specific climate model version and returns the output being calibrated to observations.
model_start_year      = First year to run the model (not necessarily first year of the calibration if model initializes earlier).
end_year              = The final year to run the model calibration (defaults to 2017).
joint_antarctic_prior = TRUE/FALSE check for whether to use a joint normal prior distribution (TRUE = option 1 described
                        above) or fitted marginal kernel density estimates (FLASE = option 2 described above).
source
MimiBRICK.construct_brick_log_priorMethod
construct_brick_log_prior(joint_antarctic_prior::Bool; calibration_data_dir::Union{String, Nothing} = nothing)

Calculate total (log) prior probability for brick.

Description: This creates a function that will calculate the total (log) prior probability of the uncertain model, initial condition, and statistical process parameters specific to the standalone BRICK model. It uses non-uniform priors for the Antarctic ice sheet parameters, informed by a previous model calibration to paleo data. There are two options for the Antarctic priors (1) fitting marginal distributions using a kernel density estimation or (2) fitting a multivariate normal distribution that accounts for correlations that emerge during the paleo calibration (note, many of the marginal paleo pdfs are not normally distributed).

Function Arguments:

  joint_antarctic_prior   = TRUE/FALSE check for whether to use a joint normal prior distribution (TRUE = option 1 described
                          above) or fitted marginal kernel density estimates (FLASE = option 2 described above).
  calibration_data_dir    = Data directory for calibration data. Defaults to package calibration data directory, 
                            changing this is not recommended.
source
MimiBRICK.construct_doeclimbrick_log_posteriorMethod
construct_doeclimbrick_log_posterior(f_run_model!; model_start_year::Int=1850, calibration_end_year::Int=2017, joint_antarctic_prior::Bool=false, uniform_ECS::Bool=false)

Calculate log posterior for doeclimbrick.

Description: This creates a function that calculates the log-posterior probability of the uncertain model, initial condition, and statistical process parameters.

Function Arguments:

f_run_model           = A function that runs the specific climate model version and returns the output being calibrated to observations.
model_start_year      = First year to run the model (not necessarily first year of the calibration if model initializes earlier).
end_year              = The final year to run the model calibration (defaults to 2017).
joint_antarctic_prior = TRUE/FALSE check for whether to use a joint normal prior distribution (TRUE = option 1 described
                        above) or fitted marginal kernel density estimates (FLASE = option 2 described above).
uniform_ECS           = TRUE/FALSE check for whether or not to use a uniform prior distribution for the equilibrium
                        climate sensitivity (true = use uniform).
source
MimiBRICK.construct_doeclimbrick_log_priorMethod
construct_doeclimbrick_log_prior(joint_antarctic_prior::Bool, uniform_ECS::Bool; calibration_data_dir::Union{String, Nothing} = nothing)

Calculate total (log) prior probability for doeclimbrick.

Description: This creates a function that will calculate the total (log) prior probability of the uncertain model, initial condition, and statistical process parameters specific to the DOECLIM-BRICK model. It uses non-uniform priors for the Antarctic ice sheet parameters, informed by a previous model calibration to paleo data. There are two options for the Antarctic priors (1) fitting marginal distributions using a kernel density estimation or (2) fitting a multivariate normal distribution that accounts for correlations that emerge during the paleo calibration (note, many of the marginal paleo pdfs are not normally distributed).

Function Arguments:

  joint_antarctic_prior   = TRUE/FALSE check for whether to use a joint normal prior distribution (TRUE = option 1 described
                          above) or fitted marginal kernel density estimates (FLASE = option 2 described above).
  uniform_ECS             = TRUE/FALSE check for whether or not to use a uniform prior distribution for the equilibrium
                          climate sensitivity (true = use uniform).
  calibration_data_dir    = Data directory for calibration data. Defaults to package calibration data directory, changing this is not recommended.
source
MimiBRICK.construct_run_brickMethod
construct_run_brick(calibration_start_year::Int, calibration_end_year::Int)

Create a function to run the BRICK model over the historic period.

source
MimiBRICK.construct_run_sneasybrickMethod
construct_run_sneasybrick(calibration_start_year::Int, calibration_end_year::Int)

Create a function to run the SNEASY+BRICK model over the historic period.

source
MimiBRICK.construct_sneasybrick_log_posteriorMethod
construct_sneasybrick_log_posterior(f_run_model!; model_start_year::Int=1850, calibration_end_year::Int=2017, joint_antarctic_prior::Bool=false, uniform_ECS::Bool=false)

Calculate log posterior for sneasybrick.

Description: This creates a function that calculates the log-posterior probability of the uncertain model, initial condition, and statistical process parameters.

Function Arguments:

f_run_model           = A function that runs the specific climate model version and returns the output being calibrated to observations.
model_start_year      = First year to run the model (not necessarily first year of the calibration if model initializes earlier).
end_year              = The final year to run the model calibration (defaults to 2017).
joint_antarctic_prior = TRUE/FALSE check for whether to use a joint normal prior distribution (TRUE = option 1 described
                        above) or fitted marginal kernel density estimates (FLASE = option 2 described above).
uniform_ECS           = TRUE/FALSE check for whether or not to use a uniform prior distribution for the equilibrium
                        climate sensitivity (true = use uniform).
source
MimiBRICK.construct_sneasybrick_log_priorMethod
construct_sneasybrick_log_prior(joint_antarctic_prior::Bool, uniform_ECS::Bool; calibration_data_dir::Union{String, Nothing} = nothing)

Calculate total (log) prior probability for sneasybrick.

Description: This creates a function that will calculate the total (log) prior probability of the uncertain model, initial condition, and statistical process parameters specific to the SNEASY-BRICK model. It uses non-uniform priors for the Antarctic ice sheet parameters, informed by a previous model calibration to paleo data. There are two options for the Antarctic priors (1) fitting marginal distributions using a kernel density estimation or (2) fitting a multivariate normal distribution that accounts for correlations that emerge during the paleo calibration (note, many of the marginal paleo pdfs are not normally distributed).

Function Arguments:

  joint_antarctic_prior   = TRUE/FALSE check for whether to use a joint normal prior distribution (TRUE = option 1 described
                          above) or fitted marginal kernel density estimates (FLASE = option 2 described above).
                          = TRUE/FALSE check for whether or not to use a uniform prior distribution for the equilibrium
                          climate sensitivity (true = use uniform).
  calibration_data_dir    = Data directory for calibration data. Defaults to package calibration data directory, changing this is not recommended.
source
MimiBRICK.create_brick_doeclimMethod
create_brick_doeclim(;rcp_scenario::String = "RCP85", start_year::Int=1850, end_year::Int=2020)

Return a Mimi model instance with MimiBRICK and DOECLIM coupled together.

Description: This function loads forcing data, sets up model parameters, and makes the model component variable connections.

Function Arguments:

    rcp_scenario = RCP scenario for exogenous forcing
    start_year   = initial year of the simulation period
    end_year     = ending year of the simulation period
source
MimiBRICK.create_sneasy_brickMethod
create_sneasy_brick(;rcp_scenario::String = "RCP85", start_year::Int=1850, end_year::Int=2020)

Return a Mimi model instance with MimiBRICK and MimiSNEASY coupled together.

Description: This function loads forcing data, sets up model parameters, and makes the model component variable connections.

Function Arguments:

    rcp_scenario = RCP scenario for exogenous forcing
    start_year   = initial year of the simulation period
    end_year     = ending year of the simulation period
source
MimiBRICK.downscale_brickMethod
downscale_brick(;lon::Float64, 
                            lat::Float64, 
                            results_dir::String, 
                            proj_or_hind::String, 
                            ensemble_or_map::String, 
                            model_config::String, 
                            rcp_scenario::String="RCP85"
                        )

Downscale BRICK projections to a single point, using either the whole ensemble or only the maximum a posteriori ensemble member. Note this function assumes a specific folder structure and file naming within the top level results_dir.

Function Arguments:

- lon = longitude (degrees East) of location for downscaling
- lat = latitude (degrees North) of location for downscaling
- results_dir = the top level directory of results
- proj_or_hind = "proj" for projections, or "hind" for hindcast
- ensemble_or_map = "ensemble" for entire posterior ensemble, or "map" for the maximum a posteriori ensemble member (single simulation)
- model_config = "brick", "doeclimbrick", or "sneasybrick"
- rcp_scenario = "RCP26", "RCP45", "RCP60", or "RCP85" (default). Doesn't matter for hindcast.
source
MimiBRICK.get_modelMethod
get_model(;rcp_scenario::String="RCP85", start_year::Int=1850, end_year::Int=2020)

Return a MimiBRICK model instance that can be modified and run.

Function Arguments:

  rcp_scenario = RCP scenario for exogenous forcing
  start_year   = initial year of the simulation period
  end_year     = ending year of the simulation period
source
MimiBRICK.hetero_logl_ar1Method
hetero_logl_ar1(residuals::Array{Float64,1}, σ::Float64, ρ::Float64, ϵ::Array{Union{Float64, Missings.Missing},1})

Calculate AR(1) log-likelihood.

Description: This function calculates the AR(1) log-likelihood in terms of the data-model residuls when accounting for time-varying observation errors. It follows "The Effects of Time-Varying Observation Errors on Semi-Empirical Sea-Level Projections" (Ruckert et al., 2017) DOI 10.1007/s10584-016-1858-z.

Function Arguments:

  residuals = A vector of data-model residuals.
  σ         = AR(1) innovation standard deviation.
  ρ         = AR(1) autocorrelation term.
  ϵ         = A vector of time-varying observation error estimates (from calibration data sets).
source
MimiBRICK.hetero_logl_car1Method
hetero_logl_car1(residuals::Array{Float64,1}, indices::Array{Int64,1}, σ²_white_noise::Float64, α₀::Float64, ϵ::Array{Union{Float64, Missings.Missing},1})

Calculate CAR(1) log-likelihood.

Description: This function calculates the continuous time autoregressive, or CAR(1), log-likelihood for irregularly spaced data in terms of the data-model residuls when accounting for time-varying observation errors. It builds off of "The Effects of Time-Varying Observation Errors on Semi-Empirical Sea-Level Projections" (Ruckert et al., 2017) DOI 10.1007/s10584-016-1858-z and "The Analysis of Irregularly Observed Stochastic Astronomical Time-Series – I. Basics of Linear Stochastic Differential Equations" (Koen, 2005) doi.org/10.1111/j.1365-2966.2005.09213.x

Function Arguments:

  residuals      = A vector of data-model residuals.
  indices        = Index positions of observations relative to model time horizon (i.e. the first model time period = 1, the second = 2, etc.).
  σ²_white_noise = Variance of the continuous white noise process.
  α₀             = Parameter describing correlation memory of CAR(1) process.
  ϵ              = A vector of time-varying observation error estimates (from calibration data sets).
source
MimiBRICK.load_calibration_dataMethod
load_calibration_data(model_start_year::Int, last_calibration_year::Int; last_sea_level_norm_year::Int=1990, calibration_data_dir::Union{Nothing, String} = nothing)

Load and clean up data used for model calibration.

Description: This function loads, cleans up, and merges all of the calibration data into a single dataframe.

Function Arguments:

- model_start_year         = The first year to include in the calibration data set.
- last_calibration_year    = The last year to run the model for calibration (i.e. 1980 will not consider any post-1980 observations).
- last_sea_level_norm_year = Some sea level data sets may need to be normalized to different years depending on when the calibration ends (this
                             may be necessary for out-of-sample tests). These data sets will be normalized from 1961-last norm year, default = 1961-1990.
- calibration_data_dir    = Data directory for calibration data. Defaults to package calibration data directory, changing this is not recommended.
source
MimiBRICK.next_latMethod
next_lat(lat::Float64, inc::Int, direction::Symbol)

Increment latitude by inc in either positive direction (direction=:increase) or in the negative direction (direction=:decrease). Assumes latitude runs from -90 to 90 (deg N).

source
MimiBRICK.next_lonMethod
next_lon(lon::Float64, inc::Int, direction::Symbol)

Increment longitude by inc in either positive direction (direction=:increase) or in the negative direction (direction=:decrease). Assumes longitude runs from 0 to 360 (deg E).

source
MimiBRICK.replicate_errorsMethod
replicate_errors(start_year::Int, end_year::Int, error_data)

REPLICATE TIME-VARYING OBSERVATION ERRORS FOR PERIODS WITHOUT COVERAGE.

Description: This function creates a time-series of observation errors for the entire model time horizon. For years without observation error estimates, the error remains constant at the average of the ten nearest error values in time.

Function Arguments:

  start_year = The first year to run the climate model.
  end_year   = The final year to run the climate model.
  error_data = A vector of time-varying observation errors supplied with each calibration data set.
source
MimiBRICK.run_calibrationMethod
run_calibration(;   output_dir::String, 
                    model_config="brick", 
                    calibration_start_year=1850, 
                    calibration_end_year=2005,
                    total_chain_length=1000, 
                    burnin_length=0, 
                    threshold_gr=1.1, 
                    num_walkers=2,
                    size_subsample=1000, 
                    start_from_priors=false,
                    calibration_data_dir::Union{String, Nothing} = nothing
                )

This function carries out a Markov chain Monte Carlo calibration of BRICK. This includes one of the following possible model configurations set with model_config (1) BRICK standalone (forced by input global mean surface temperatures and ocean heat uptake data) (2) DOECLIM+BRICK (3) SNEASY+BRICK

source
MimiBRICK.run_hindcastMethod
run_hindcast(; output_dir::String,
                model_config::String = "brick",
                start_year::Int = 1850,
                end_year = 2017,
            )

Function to run BRICK (standalone, or with DOECLIM or SNEASY) over the historic period and save the hindcast results to CSV files.

Function Arguments:

- outdir - paths for results files - subsample of model parameters, and associated log-posterior scores, and printed results of this function
- model_config (default = "brick") - model configuration with possible options: (1) "brick", (2) "doeclimbrick", (3) "sneasybrick"
- start_year (default = 1850) - start year for calibration
- end_year (default = 2017) - end year for calibration
source
MimiBRICK.run_projectionsMethod
run_projections(; output_dir::String,
                    model_config::String = "brick",
                    rcp_scenario::String = "RCP85",
                    start_year::Int = 1850,
                    end_year = 2300,
                )

Function to run BRICK (standalone, or with DOECLIM or SNEASY) over the projections period and save the results to CSV files.

Function Arguments:

- outdir - paths for results files - subsample of model parameters, and associated log-posterior scores, and printed results of this function
- model_config (default = "brick") - model configuration with possible options: (1) "brick", (2) "doeclimbrick", (3) "sneasybrick"
- rcp_scenario (default = "RCP85) - RCP scenario with possible options: (1) RCP26, (2) RCP45, (3) RCP60, (4) RCP85
- start_year (default = 1850) - start year for calibration
- end_year (default = 2300) - end year for calibration
source
MimiBRICK.simulate_ar1_noiseMethod
simulate_ar1_noise(n::Int, σ::Float64, ρ::Float64, ϵ::Array{Float64,1})

Simulate stationary AR(1) process with time varying observation errors.

Description: This function simulates a stationary AR(1) process (given time-varying observation errors supplied with each calibration data set) to superimpose noise onto the climate model projections.

Function Arguments:

  n = Number of time periods (years) the model is being run for.
  σ = Calibrated standard deviation.
  ρ = Calibrated autocorrelation coefficient.
  ϵ = Time-varying observation errors.
source
MimiBRICK.simulate_car1_noiseMethod
simulate_car1_noise(n, α₀, σ²_white_noise, ϵ)

Simulate stationary CAR(1) process with time varying observation errors.

Description: This function simulates a stationary CAR(1) process (given time-varying observation errors supplied with each calibration data set) to superimpose noise onto the climate model projections.

Function Arguments:

  n              = Number of time periods (years) the model is being run for.
  α₀             = Calibrated term describing correlation memory of CAR(1) process.
  σ²_white_noise = Calibrated continuous white noise process variance term.
  ϵ              = Time-varying observation errors.
source
MimiBRICK.truncated_kernelMethod
truncated_kernel(data, lower_bound, upper_bound)

Calculate kernel density estimates with truncated bounds.

Description: This function creates a fitted kernel density object (from KernelDensity.jl package) and crops the edges to user-specified lower and upper bounds (setting boundaries beforehand, can lead to issues with some parameters because, "Due to the fourier transforms used internally, there should be sufficient spacing to prevent wrap-around at the boundaries")

Function Arguments:

  data        = Parameter samples used for the kernel density estimate.
  lower_bound = Minimum parameter value corresponding to the lower bound of the fitted density.
  upper_bound = Maximum parameter value corresponding to the upper bound of the fitted density.
source

Author Contributions

  • TW: initial model development, software development, model calibration and validation, conceptualization, projection direction and overall management
  • LR: software development, package maintenance, conceptualization
  • FE: software development, model calibration and validation, conceptualization
  • VS: software testing, model calibration and validation, conceptualization
  • AB: initial model development, software testing, conceptualization
  • KK: software testing, conceptualization
  • DA: software development, package maintenance, conceptualization, project direction

License

Copyright 2022 Tony Wong, Lisa Rennels, Frank Errickson, Vivek Srikrishnan, Alexander Bakker, Klaus Keller, and David Anthoff

This file and codes in this repository are part of MimiBRICK.jl (Building blocks for Relevant Ice and Climate Knowledge). MimiBRICK.jl is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

MimiBRICK.jl is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with MimiBRICK.jl (LICENSE.md)). If not, see http://www.gnu.org/licenses/.