Parameters

There are different types of parameters in ODINN, holding specific information for different modelling aspects. All the types of parameters are wrapped into a parent Parameter type, which is threaded throughout ODINN.jl.

Sleipnir.ParametersType
mutable struct Parameters{PPHY <: AbstractEmptyParams, PSIM <: AbstractEmptyParams, PHY <: AbstractEmptyParams,
    PSOL <: AbstractEmptyParams, PUDE <: AbstractEmptyParams, PINV <: AbstractEmptyParams}

A mutable struct that holds various parameter sets for different aspects of a simulation or model.

Fields

  • physical::PPHY: Physical parameters.
  • simulation::PSIM: Simulation parameters.
  • hyper::PHY: Hyperparameters.
  • solver::PSOL: Solver parameters.
  • UDE::PUDE: Universal Differential Equation (UDE) parameters.
  • inversion::PINV: Inversion parameters.

Type Parameters

  • PPHY: Type of the physical parameters, must be a subtype of AbstractEmptyParams.
  • PSIM: Type of the simulation parameters, must be a subtype of AbstractEmptyParams.
  • PHY: Type of the hyperparameters, must be a subtype of AbstractEmptyParams.
  • PSOL: Type of the solver parameters, must be a subtype of AbstractEmptyParams.
  • PUDE: Type of the UDE parameters, must be a subtype of AbstractEmptyParams.
  • PINV: Type of the inversion parameters, must be a subtype of AbstractEmptyParams.
source
ODINN.ParametersFunction

Constructor for the Parameters type. Since some of the subtypes of parameters are defined in different packages of the ODINN ecosystem, this constructor will call the constructors of the different subtypes and return a Parameters object with the corresponding subtypes. The Parameters mutable struct is defined in Sleipnir.jl using abstract types, which are later on defined in the different packages of the ODINN ecosystem.

Parameters(;
        physical::PhysicalParameters = PhysicalParameters(),
        simulation::SimulationParameters = SimulationParameters(),
        solver::SolverParameters = SolverParameters(),
        hyper::Hyperparameters = Hyperparameters(),
        UDE::UDEparameters = UDEparameters()
        inversion::InversionParameters = InversionParameters()
        )

Keyword arguments

  • physical::PhysicalParameters: Physical parameters for the simulation.
  • simulation::SimulationParameters: Parameters related to the simulation setup.
  • solver::SolverParameters: Parameters for the solver configuration.
  • hyper::Hyperparameters: Hyperparameters for the model.
  • UDE::UDEparameters: Parameters specific to the UDE (Universal Differential Equation).
  • inversion::InversionParameters: Parameters for inversion processes.
source

The main child types of parameters are the following ones:

Simulation parameters

Simulation parameters are used to specify anything related to ODINN simulations, ranging from types, working directories to multiprocessing.

Sleipnir.SimulationParametersType

A structure to hold simulation parameters for a simulation in ODINN.

struct SimulationParameters{I <: Integer, F <: AbstractFloat, VM <: VelocityMapping} <: AbstractParameters

Fields

  • use_MB::Bool: Flag to indicate whether mass balance should be used.
  • use_iceflow::Bool: Flag to indicate whether ice flow should be used.
  • plots::Bool: Flag to indicate whether plots should be generated.
  • use_velocities::Bool: Flag to indicate whether velocities should be calculated.
  • overwrite_climate::Bool: Flag to indicate whether to overwrite climate data.
  • use_glathida_data::Bool: Flag to indicate whether to use GLATHIDA data.
  • tspan::Tuple{F, F}: Time span for the simulation.
  • step::F: Time step for the simulation.
  • multiprocessing::Bool: Flag to indicate whether multiprocessing should be used.
  • workers::I: Number of workers for multiprocessing.
  • working_dir::String: Directory for working files.
  • test_mode::Bool: Flag to indicate whether to run in test mode.
  • rgi_paths::Dict{String, String}: Dictionary of RGI paths.
  • ice_thickness_source::String: Source of ice thickness data.
  • mapping::VM: Mapping to use in order to grid the data from the coordinates of the velocity product datacube to the glacier grid.
  • gridScalingFactor::I: Grid downscaling factor, used to speed-up the tests. Default value is 1 which means no downscaling is applied.
source
Sleipnir.SimulationParametersMethod

Constructor for SimulationParameters type, including default values.

SimulationParameters(;
    use_MB::Bool = true,
    use_iceflow::Bool = true,
    plots::Bool = true,
    use_velocities::Bool = true,
    overwrite_climate::Bool = false,
    use_glathida_data::Bool = false,
    tspan::Tuple{F, F} = (2010.0,2015.0),
    step::F = 1/12,
    multiprocessing::Bool = true,
    workers::I = 4,
    working_dir::String = "",
    test_mode::Bool = false,
    rgi_paths::Dict{String, String} = Dict{String, String}(),
    ice_thickness_source::String = "Farinotti19",
    mapping::VM = MeanDateVelocityMapping(),
    gridScalingFactor::I = 1,
) where {I <: Integer, F <: AbstractFloat, VM <: VelocityMapping}

Keyword arguments

  • use_MB::Bool: Whether to use mass balance (default: true).
  • use_iceflow::Bool: Whether to use ice flow (default: true).
  • plots::Bool: Whether to generate plots (default: true).
  • use_velocities::Bool: Whether to calculate velocities (default: true).
  • overwrite_climate::Bool: Whether to overwrite climate data (default: false).
  • use_glathida_data::Bool: Whether to use GLATHIDA data (default: false).
  • float_type::DataType: Data type for floating point numbers (default: Float64).
  • int_type::DataType: Data type for integers (default: Int64).
  • tspan::Tuple{F, F}: Time span for the simulation (default: (2010.0, 2015.0)).
  • step::F: Time step for the simulation (default: 1/12).
  • multiprocessing::Bool: Whether to use multiprocessing (default: true).
  • workers::I: Number of workers for multiprocessing (default: 4).
  • working_dir::String: Working directory for the simulation (default: "").
  • test_mode::Bool: Whether to run in test mode (default: false).
  • rgi_paths::Dict{String, String}: Dictionary of RGI paths (default: Dict{String, String}()).
  • ice_thickness_source::String: Source of ice thickness data, either "Millan22" or "Farinotti19" (default: "Farinotti19").
  • mapping::VM: Mapping to use in order to grid the data from the coordinates of the velocity product datacube to the glacier grid.
  • gridScalingFactor::I: Grid downscaling factor, used to speed-up the tests. Default value is 1 which means no downscaling is applied.

Returns

  • simulation_parameters: A new SimulationParameters object.

Throws

  • AssertionError: If ice_thickness_source is not "Millan22" or "Farinotti19".

Notes

  • If the global variable ODINN_OVERWRITE_MULTI is set to true, multiprocessing is is enabled in any case and the number of workers specified in the simulation parameters must correspond to the number of processes with which Julia has been started. This is to allow the documentation to build successfully in ODINN as we cannot change the number of process in the CI.
source

Physical parameters

Physical parameters are used to store physical constants and variables used in the physical and machine learning models.

Sleipnir.PhysicalParametersType

A structure representing physical parameters used in simulations.

PhysicalParameters{F <: AbstractFloat}

Fields

  • ρ::F: Density of ice.
  • g::F: Gravitational acceleration.
  • ϵ::F: Regularization used in the square root of norms for AD numerical stability.
  • η₀::F: Initial viscosity.
  • maxA::F: Maximum value for A (Glen's coefficient).
  • minA::F: Minimum value for A (Glen's coefficient).
  • maxTlaw::F: Maximum value of Temperature used in simulations on fake law.
  • minTlaw::F: Minimum value of Temperature used in simulations on fake law.
  • noise_A_magnitude::F: Magnitude of noise in A.
source
Sleipnir.PhysicalParametersMethod

Initialize the physical parameters of a model.

PhysicalParameters(;
    ρ::Float64 = 900.0,
    g::Float64 = 9.81,
    ϵ::Float64 = 1e-10,
    η₀::F = 1.0,
    maxA::Float64 = 8e-17,
    minA::Float64 = 8.5e-20,
    maxTlaw::Float64 = 1.0,
    minTlaw::Float64 = -25.0,
    noise_A_magnitude::Float64 = 5e-18
    )

Keyword arguments

- `ρ`: Density of ice.
- `g`: Gravitational acceleration.
- `ϵ`: Regularization used in the square root of norms for AD numerical stability.
- `η₀`: Initial viscosity.
- `maxA`: Maximum value for `A` (Glen's coefficient).
- `minA`: Minimum value for `A` (Glen's coefficient).
- `maxTlaw`: Maximum value of Temperature used in simulations on fake law.
- `minTlaw`: Minimum value of Temperature used in simulations on fake law.
- `noise_A_magnitude`: Magnitude of noise added to A
source

Solver parameters

Solver parameters determine all aspects related to the numerical scheme used to solve the differential equations of glacier ice flow.

Huginn.SolverParametersType

A mutable struct that holds parameters for the solver.

SolverParameters{F <: AbstractFloat, I <: Integer}

Fields

  • solver::OrdinaryDiffEq.OrdinaryDiffEqAdaptiveAlgorithm: The algorithm used for solving differential equations.
  • reltol::F: The relative tolerance for the solver.
  • step::F: The step size for the solver.
  • tstops::Union{Nothing, Vector{F}}: Optional vector of time points where the solver should stop for the callbacks.
  • save_everystep::Bool: Flag indicating whether to save the solution at every step.
  • progress::Bool: Flag indicating whether to show progress during the solving process.
  • progress_steps::I: The number of steps between progress updates.
  • maxiters::I: Maximum number of iterations to perform in the iceflow solver.
source
Huginn.SolverParametersMethod

Constructs a SolverParameters object with the specified parameters or using default values.

SolverParameters(; solver::OrdinaryDiffEq.OrdinaryDiffEqAdaptiveAlgorithm = RDPK3Sp35(),
                  reltol::F = 1e-12,
                  step::F = 1.0/12.0,
                  tstops::Union{Nothing,Vector{F}} = nothing,
                  save_everystep = false,
                  progress::Bool = true,
                  progress_steps::I = 10,
                  maxiters::I = Int(1e5),
                ) where {F <: AbstractFloat, I <: Integer}

Arguments

  • solver::OrdinaryDiffEq.OrdinaryDiffEqAdaptiveAlgorithm: The ODE solver algorithm to use. Defaults to RDPK3Sp35().
  • reltol::F: The relative tolerance for the solver. Defaults to 1e-12.
  • step::F: The step size for the callbacks. These are mainly used to run the surface mass balance model. Defaults to 1.0/12.0 (i.e. a month).
  • tstops::Union{Nothing, Vector{F}}: Optional vector of time points where the solver should stop. Defaults to nothing.
  • save_everystep::Bool: Whether to save the solution at every step. Defaults to false.
  • progress::Bool: Whether to show progress during the solving process. Defaults to true.
  • progress_steps::I: The number of steps between progress updates. Defaults to 10.
  • maxiters::I: Maximum number of iterations to perform in the iceflow solver. Defaults to 1e5.

Returns

  • solver_parameters: A SolverParameters object constructed with the specified parameters.
source

Hyperparameters

Hyperparameters determine different aspects of a given machine learning model. For now, these are focused on neural networks, but we plan to extend them in the future for other types of regressors.

ODINN.HyperparametersType
mutable struct Hyperparameters{F <: AbstractFloat, I <: Integer} <: AbstractParameters

A mutable struct that holds hyperparameters for training a machine learning model.

Keyword arguments

  • current_epoch::I: The current epoch number.
  • current_minibatch::I: The current minibatch number.
  • loss_history::Vector{F}: A vector storing the history of loss values.
  • optimizer::Union{Optim.FirstOrderOptimizer, Flux.Optimise.AbstractOptimiser, Optimisers.AbstractRule}: The optimizer used for training.
  • loss_epoch::F: The loss value for the current epoch.
  • epochs::I: The total number of epochs for training.
  • batch_size::I: The size of each minibatch.
source
ODINN.HyperparametersMethod
Hyperparameters(; current_epoch::Int64 = 1, current_minibatch::Int64 = 1, loss_history::Vector{Float64} = Vector{Float64}(), optimizer::Union{Optim.FirstOrderOptimizer, Flux.Optimise.AbstractOptimiser, Optimisers.AbstractRule} = BFGS(initial_stepnorm=0.001), loss_epoch::Float64 = 0.0, epochs::Int64 = 50, batch_size::Int64 = 15)

Constructs a Hyperparameters object with the specified parameters.

Arguments

  • current_epoch::Int64: The current epoch number. Defaults to 1.
  • current_minibatch::Int64: The current minibatch number. Defaults to 1.
  • loss_history::Vector{Float64}: A vector to store the history of loss values. Defaults to an empty vector.
  • optimizer::Union{Optim.FirstOrderOptimizer, Flux.Optimise.AbstractOptimiser, Optimisers.AbstractRule}: The optimizer to be used. Defaults to BFGS(initial_stepnorm=0.001).
  • loss_epoch::Float64: The loss value for the current epoch. Defaults to 0.0.
  • epochs::Int64: The total number of epochs. Defaults to 50.
  • batch_size::Int64: The size of each minibatch. Defaults to 15.

Returns

  • A Hyperparameters object initialized with the provided values.
source

UDE parameters

Universal Differential Equation (UDE) parameters are used to determine different modelling choices regarding the use of UDEs, such as wich sensitivity algorithm or optimization method to use.

ODINN.UDEparametersType

A mutable struct that holds parameters for a UDE (Universal Differential Equation).

UDEparameters{ADJ <: AbstractAdjointMethod} <: AbstractParameters

Fields

  • sensealg::SciMLBase.AbstractAdjointSensitivityAlgorithm: The sensitivity algorithm used for adjoint sensitivity analysis.
  • optimization_method::String: The optimization method to be used.
  • target::Symbol: The target variable for the optimization.
source
ODINN.UDEparametersMethod
UDEparameters(; sensealg, optim_autoAD, grad, optimization_method, empirical_loss_function, target) where {ADJ <: AbstractAdjointMethod}

Create a UDEparameters object for configuring the sensitivity analysis and optimization of a Universal Differential Equation (UDE).

Keyword Arguments

  • sensealg::SciMLBase.AbstractAdjointSensitivityAlgorithm: The sensitivity algorithm to use for adjoint calculations. Defaults to GaussAdjoint(autojacvec=SciMLSensitivity.EnzymeVJP()).
  • optim_autoAD::AbstractADType: The automatic differentiation type for optimization. Defaults to Optimization.AutoEnzyme().
  • grad::ADJ: The adjoint gradient computation method. Defaults to SciMLSensitivityAdjoint().
  • optimization_method::String: The optimization method to use. Must be either "AD+AD" or "AD+Diff". Defaults to "AD+AD".
  • empirical_loss_function::AbstractLoss: The loss function to use for optimization. Defaults to LossH().
  • target::Union{Symbol, Nothing}: The target variable for optimization. Defaults to :A.

Returns

  • A UDEparameters object configured with the specified sensitivity, optimization, and loss settings.

Description

This function creates a UDEparameters object that encapsulates the configuration for sensitivity analysis, optimization, and loss computation in a Universal Differential Equation (UDE) framework. It verifies that the provided optimization_method is valid and constructs the solver parameters accordingly.

Notes

  • The optimization_method must be either "AD+AD" (automatic differentiation for both forward and backward passes) or "AD+Diff" (automatic differentiation combined with finite differences).
  • The empirical_loss_function determines how the loss is computed during optimization.
source