API
UniversalDiffEq.BayesianUDE
— TypeBayesianUDE
Basic data structure used to the model structure, parameters and data for Bayesian UDE and NODE models. ...
Elements
- times: a vector of times for each observation
- data: a matrix of observations at each time point
- X: a DataFrame with any covariates used by the model
- data_frame: a DataFrame with columns for the time of each observation and values of the state variables
- parameters: a ComponentArray that stores model parameters
- loss_function: the loss function used to fit the model
- process_model: a Julia mutable struct used to define model predictions
- process_loss: a Julia mutable struct used to measure the performance of model predictions
- observation_model: a Julia mutable struct used to predict observations given state variable estimates
- observation_loss: a Julia mutable struct used to measure the performance of the observation model
- process_regularization: a Julia mutable struct used to store data needed for process model regularization
- observation_regularization: a Julia mutable struct used to store data needed for observation model regularization
- constructor: A function that initializes a UDE model with identical structure.
...
UniversalDiffEq.LossFunction
— TypeLossFunction
A Julia mutable struct that stores the loss function and parameters. ...
Elements
- parameters: ComponentArray
- loss: Function
...
UniversalDiffEq.MultiProcessModel
— TypeMultiProcessModel
A Julia mutable struct that stores the functions and parameters for the process model. ...
Elements
parameters
: ComponentArraypredict
: Function that predicts one time step aheadforecast
: Function that is a modified version of predict to improve performace when extrapolatingcovariates
: Function that returns the values of the covariates at each point in timeright_hand_side
: Function that returns the right-hand side of a differential equation (i.e., the relationships between state variables and parameters)
...
UniversalDiffEq.MultiTimeSeriesNetwork
— MethodMultiTimeSeriesNetwork(inputs,series,outputs; kwargs...)
Builds a neural network object and returns randomly initialized parameter values. The neural network object can be evaluated as a function that takes three arguments: a vector x
, an integer i
, and a named tuple with the neural network weights and biases parameters
.
# kwargs
hidden
: the number of neurons in the hidden layer. The default is 10.nonlinearity
: the activation function used in the neural network. The default is the hyperbolic tangent functiontanh
.seed
: random number generator seed for initializing the neural network weights.distance
: determines the level of difference between the functions estimated for each time series. The default is 1.0.
UniversalDiffEq.ProcessModel
— TypeProcessModel
A Julia mutable struct that stores the functions and parameters for the process model. ...
Elements
- parameters: ComponentArray
- predict: Function the predict one time step ahead
- forecast: Function, a modified version of predict to improve performance when extrapolating
- covariates: Function that returns the value of the covariates at each point in time.
...
UniversalDiffEq.Regularization
— TypeRegularization
A Julia mutable struct that stores the loss function and parameters. ...
Elements
- reg_parameters: ComponentArray
- loss: Function
...
UniversalDiffEq.SimpleNeuralNetwork
— MethodSimpleNeuralNetwork(inputs,outputs; kwargs ...)
Builds a neural network object and returns randomly initialized parameter values. The neural network object can be evaluated as a function that takes two arguments: a vector x
and a named tuple with the neural network weights and biases parameters
.
# kwargs
hidden
: the number of neurons in the hidden layer. The default is 10.nonlinearity
: the activation function used in the neural network. The default is the hyperbolic tangent functiontanh
.seed
: random number generator seed for initializing the neural network weights.
UniversalDiffEq.UDE
— TypeUDE
Basic data structure used to the model structure, parameters and data for UDE and NODE models. ...
Elements
- times: a vector of times for each observation
- data: a matrix of observations at each time point
- X: a DataFrame with any covariates used by the model
- data_frame: a DataFrame with columns for the time of each observation and values of the state variables
- parameters: a ComponentArray that stores model parameters
- loss_function: the loss function used to fit the model
- process_model: a Julia mutable struct used to define model predictions
- process_loss: a Julia mutable struct used to measure the performance of model predictions
- observation_model: a Julia mutable struct used to predict observations given state variable estimates
- observation_loss: a Julia mutable struct used to measure the performance of the observation model
- process_regularization: a Julia mutable struct used to store data needed for process model regularization
- observation_regularization: a Julia mutable struct used to store data needed for observation model regularization
- constructor: A function that initializes a UDE model with identical structure.
- timecolumnname: A string with the name of the column used for
- weights
- variablecolumnname
- valuecolumnname
...
UniversalDiffEq.BFGS!
— Method BFGS!(UDE, kwargs ...)
Minimizes the loss function of the UDE
model using the BFGS algorithm is the initial step norm equal to initial_step_norm
. The function will print the value of the loss function after each iteration when verbose
is true.
kwargs
initial_step_norm
: Initial step norm for BFGS algorithm. Default is0.01
.verbose
: Should the training loss values be printed? Default isfalse
.
UniversalDiffEq.BayesianCustomDerivatives
— MethodBayesianCustomDerivatives(data::DataFrame,X,derivs!::Function,initial_parameters;kwargs ... )
When a data frame X
is supplied the model will run with covariates. The argument X
should have a column for time t
with the value for time in the remaining columns. The values in X
will be interpolated with a linear spline for value of time not included in the data frame.
When X
is provided the derivs function must have the form derivs!(du,u,x,p,t)
where x
is a vector with the value of the covariates at time t
.
kwargs
time_column_name
: Name of column indata
andX
that corresponds to time. Default is"time"
.variable_column_name
: Name of column inX
that corresponds to the variables. Default isnothing
.value_column_name
: Name of column inX
that corresponds to the covariates. Default isnothing
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.
UniversalDiffEq.BayesianCustomDerivatives
— MethodBayesianCustomDerivatives(data::DataFrame,derivs!::Function,initial_parameters;kwargs ... )
Constructs a Bayesian UDE model for the data set data
based on user defined derivatives derivs
. An initial guess of model parameters are supplied with the initial_parameters
argument.
...
Arguments
- data: a DataFrame object with the time of observations in a column labeled
t
and the remaining columns the value of the state variables at each time point. - derivs: a Function of the form
derivs!(du,u,p,t)
whereu
is the value of the state variables,p
are the model parameters,t
is time, and du is updated with the value of the derivatives - init_parameters: A
NamedTuple
with the model parameters. Neural network parameters must be listed under the keyNN
.
kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.
...
UniversalDiffEq.BayesianNODE
— MethodBayesianNODE(data,X;kwargs ... )
When a data frame X
is supplied the model will run with covariates. The argument X
should have a column for time t
with the value for time in the remaining columns. The values in X
will be interpolated with a linear spline for values of time not included in the data frame.
kwargs
time_column_name
: Name of column indata
andX
that corresponds to time. Default is"time"
.variable_column_name
: Name of column inX
that corresponds to the variables. Default isnothing
.value_column_name
: Name of column inX
that corresponds to the covariates. Default isnothing
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.
UniversalDiffEq.BayesianNODE
— MethodBayesianNODE(data;kwargs ... )
Constructs a Bayesian continuous-time model for the data set data
using a single layer neural network to represent the system's dynamics.
kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.
UniversalDiffEq.CustomDerivatives
— MethodCustomDerivatives(data,derivs!,initial_parameters;kwargs ... )
Constructs a UDE model for the data set data
based on user-defined derivatives derivs
. An initial guess of model parameters are supplied with the initial_parameters
argument.
...
Arguments
- data: a DataFrame object with the time of observations in a column labeled
t
and the remaining columns the value of the state variables at each time point. - derivs: a Function of the form
derivs!(du,u,p,t)
whereu
is the value of the state variables,p
are the model parameters,t
is time, and du is updated with the value of the derivatives - init_parameters: A
NamedTuple
with the model parameters. Neural network parameters must be listed under the keyNN
.
kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.ode_solver
: method to aproximate solutions to the differntail equation. Defaul isTsit5()
ad_method
:method to evalaute derivatives of the ODE solver. Default isForwardDiffSensitivity()
...
UniversalDiffEq.CustomDerivatives
— MethodCustomDerivatives(data::DataFrame,X::DataFrame,derivs!::Function,initial_parameters;kwargs ... )
When a data frame X
is supplied the model will run with covariates. The argument X
should have a column for time t
with the value for time in the remaining columns. The values in X
will be interpolated with a linear spline for values of time not included in the data frame.
When X
is provided the derivs function must have the form derivs!(du,u,x,p,t)
where x
is a vector with the value of the covariates at time t
.
kwargs
time_column_name
: Name of column indata
andX
that corresponds to time. Default is"time"
.variable_column_name
: Name of column inX
that corresponds to the variables. Default isnothing
.value_column_name
: Name of column inX
that corresponds to the covariates. Default isnothing
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.ode_solver
: method to aproximate solutions to the differntail equation. Defaul isTsit5()
.ad_method
:method to evalaute derivatives of the ODE solver. Default isForwardDiffSensitivity()
.
UniversalDiffEq.CustomDerivatives
— MethodCustomDerivatives(data::DataFrame,derivs!::Function,initial_parameters,priors::Function;kwargs ... )
When a function's priors is supplied its value will be added to the loss function as a penalty term for user specified parameters. It should take the a single NamedTuple p
as an argument penalties for each paramter should be calculated by accessing p
with the period operator.
The prior function can be used to nudge the fitted model toward prior expectations for a parameter value. For example, the following function increases the loss when a parameter p.r
has a value other than 1.5, nad a second parameter p.beta
is greater than zeros.
function priors(p)
l = 0.01*(p.r - 1.5)^2
l += 0.01*(p.beta)^2
return l
end
kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.ode_solver
: method to aproximate solutions to the differntail equation. Defaul isTsit5()
ad_method
:method to evalaute derivatives of the ODE solver. Default isForwardDiffSensitivity()
UniversalDiffEq.CustomDifference
— MethodCustomDifference(data,step,initial_parameters;kwrags...)
Constructs a UDE model for the data set data
based on user defined difference equation step
. An initial guess of model parameters are supplied with the initial_parameters argument. ...
Arguments
- data: a DataFrame object with the time of observations in a column labeled
t
and the remaining columns the value of the state variables at each time point. - step: a Function of the form
step(u,t,p)
whereu
is the value of the state variables,p
are the model parameters. - init_parameters: A
NamedTuple
with the model parameters. Neural network parameters must be listed under the keyNN
.
kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.
...
UniversalDiffEq.CustomDifference
— MethodCustomDifference(data::DataFrame,X,step,initial_parameters;kwargs ... )
When a data frame X
is supplied the model will run with covariates. The argument X
should have a column for time t
with the value for time in the remaining columns. The values in X
will be interpolated with a linear spline for value of time not included in the data frame.
When X
is provided the step function must have the form step(u,x,t,p)
where x
is a vector with the value of the covariates at time t
.
# kwargs
time_column_name
: Name of column indata
andX
that corresponds to time. Default is"time"
.variable_column_name
: Name of column inX
that corresponds to the variables. Default isnothing
.value_column_name
: Name of column inX
that corresponds to the covariates. Default isnothing
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.
UniversalDiffEq.CustomDifference
— MethodCustomDifference(data::DataFrame,step,initial_parameters,priors::Function;kwargs ... )
When a function's priors is supplied its value will be added to the loss function as a penalty term for user-specified parameters. It should take the a single NamedTuple p
as an argument penalties for each parameter should be calcualted by accessing p
with the period operator.
function priors(p)
l = 0.01*(p.r - 1.5)^2
l += 0.01*(p.beta)^2
return l
end
kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.
UniversalDiffEq.CustomModel
— MethodCustomModel(data::DataFrame,X::DataFrame, derivs!::Function, initial_parameters; kwargs ...)
Constructs a UDE model from a DataFrame data
a DataFrame with covariates X
a function derivs
and an initial guess of the paramter values inital_parameters
. The modle strucutre can be further modified by the key word arguments to specify the relationship between the state variables and observations and the loss function. These areguments are discussed individuals below.
kwargs
link - A function that takes the value of the state variable u
and paramters p
and retuns and estiamte of the obseration y
linkparams - parameters for the link function, can be an empty NamedTuple if no paramters are used observationloss - loss function that describes the distance betwen the observed and estimated states observationparams - parameters for the obseraiton loss function - can be an empty named tuple of no paramters are needed processloss - loss funciton that describes the distance betwen the observed and predicted state tranistions. processlossparams - parameters for the process loss. statevariabletransform - a function that maps from the variables used in the optimizer to states variables used by the observaiton and prediction funitons. logpriors - prior probabilities for the model paramters + nerual network regularization timecolumnname - column that indexes time in the data frames valuecolumnname - the column that indicates the variabe in long formate covariates data sets variablecolumnname = the column that indicates the value of the variables in long formate covariates data sets regweight - weight given to regualrizing the neural network reg_type - funcrional form of regualrization "L1" or "L2"
UniversalDiffEq.EasyNODE
— MethodEasyNODE(data,X;kwargs ... )
When a data frame X
is supplied the model will run with covariates. The argument X
should have a column for time t
with the value for time in the remaining columns. The values in X
will be interpolated with a linear spline for values of time not included in the data frame.
kwargs
time_column_name
: Name of column indata
andX
that corresponds to time. Default is"time"
.variable_column_name
: Name of column inX
that corresponds to the variables. Default isnothing
.value_column_name
: Name of column inX
that corresponds to the covariates. Default isnothing
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.step_size
: Step size for ADAM optimizer. Default is0.05
.maxiter
: Maximum number of iterations in gradient descent algorithm. Default is500
.verbose
: Should the training loss values be printed?. Default isfalse
.
UniversalDiffEq.EasyNODE
— MethodEasyNODE(data;kwargs ... )
Constructs a pretrained continuous-time model for the data set data
using a single layer neural network to represent the system's dynamics.
kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.step_size
: Step size for ADAM optimizer. Default is0.05
.maxiter
: Maximum number of iterations in gradient descent algorithm. Default is500
.verbose
: Should the training loss values be printed?. Default isfalse
.
UniversalDiffEq.EasyUDE
— MethodEasyUDE(data,derivs!,initial_parameters;kwargs ... )
Constructs a pretrained UDE model for the data set data
based on user defined derivatives derivs
. An initial guess of model parameters are supplied with the initial_parameters
argument.
kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.step_size
: Step size for ADAM optimizer. Default is0.05
.maxiter
: Maximum number of iterations in gradient descent algorithm. Default is500
.verbose
: Should the training loss values be printed?. Default isfalse
.
UniversalDiffEq.EasyUDE
— MethodEasyUDE(data::DataFrame,X,derivs!::Function,initial_parameters;kwargs ... )
When a data frame X
is supplied the model will run with covariates. The argument X
should have a column for time t
with the value for time in the remaining columns. The values in X
will be interpolated with a linear spline for value of time not included in the data frame. When X
is provided the derivs function must have the form derivs!(du,u,x,p,t)
where x
is a vector with the value of the covariates at time t
.
# kwargs
time_column_name
: Name of column indata
andX
that corresponds to time. Default is"time"
.variable_column_name
: Name of column inX
that corresponds to the variables. Default is"variable"
.value_column_name
: Name of column inX
that corresponds to the covariates. Default is"value"
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.step_size
: Step size for ADAM optimizer. Default is0.05
.maxiter
: Maximum number of iterations in gradient descent algorithm. Default is500
.verbose
: Should the training loss values be printed?. Default isfalse
.
UniversalDiffEq.LorenzLotkaVolterra
— MethodLorenzLotkaVolterra(;kwargs)
Create a sample dataset using the Lorenz Lotka-Volterra model as its process model:
```math
rac{dx}{dt} = rx(1-rac{x}{K}) - lpha xy + gz\
rac{dy}{dt} = hetalpha xy - my\
rac{dz}{dt} = l(w-z)\
rac{dw}{dt} = z(
ho-s) 0 w\
rac{ds}{dt} = zw-eta s
```
and an observation error following a normal distribution with mean 0 and standard deviation σ_{obs}.
# kwargs
- `plot`: Does the function return a plot? Default is `true`.
- `seed`: Seed for observation error to create repeatable examples. Default is `123`.
- `datasize`: Number of time steps generated. Default is `60`.
- `T`: Maximum timespan. Default is `3.0`.
- `sigma`: Standard deviation of observation error. Default is `0.075`.
UniversalDiffEq.LotkaVolterra
— MethodLotkaVolterra(;kwargs)
Create a sample dataset using the Lotka-Volterra predator-prey model as its process model:
```math
rac{dN}{dt} = rN - lpha NP \
rac{dP}{dt} = hetalpha NP - mP
```
and an observation error following a normal distribution with mean 0 and standard deviation σ.
# kwargs
- `plot`: Does the function return a plot? Default is `true`.
- `seed`: Seed for observation error to create repeatable examples. Default is `123`.
- `datasize`: Number of time steps generated. Default is `60`.
- `T`: Maximum timespan. Default is `3.0`.
- `sigma`: Standard deviation of observation error. Default is `0.075`.
UniversalDiffEq.MultiCustomDerivatives
— MethodMultiCustomDerivatives(data,derivs!,initial_parameters;kwargs...)
Builds a UDE model that can be trianed on multiple time series simultaniously. The user defined derivatives functions must allow for an extra argument i
that indexes over the time series in the data set (e.g. derivs!(du,u,i,)
). data
is a DataFrame object with time arguments placed in a column labeled t
and a second column with a unique index for each time series. The remaining columns have observations of the state variables at each point in time and for each time series.
kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.series_column_name
: Name of column indata
that corresponds to series. Default is"series"
.variable_column_name
: Name of column indata
that corresponds to the variables. Default is"variable"
.value_column_name
: Name of column indata
that corresponds to the covariates. Default is"value"
.proc_weight
: Weight of process erroromega_{proc}
. Default is1.0
.obs_weight
: Weight of observation erroromega_{obs}
. Default is1.0
.reg_weight
: Weight of regularization erroromega_{reg}
. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.ode_solver
: method to aproximate solutions to the differntail equation. Defaul isTsit5()
.ad_method
:method to evalaute derivatives of the ODE solver. Default isForwardDiffSensitivity()
.
UniversalDiffEq.MultiNODE
— MethodMultiNODE(data,X;kwargs...)
When a dataframe X
is supplied the model will run with covariates. the argument X
should have a column for time t
with the value for time in the remaining columns. The values in X
will be interpolated with a linear spline for values of time not included in the data frame.
kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.series_column_name
: Name of column indata
that corresponds to time. Default is"series"
.variable_column_name
: Name of column indata
that corresponds to the covariates. Default is"variable"
.value_column_name
: Name of column indata
that corresponds to the covariates. Default is"value"
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process erroromega_{proc}
. Default is1.0
.obs_weight
: Weight of observation erroromega_{obs}
. Default is1.0
.reg_weight
: Weight of regularization erroromega_{reg}
. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.ode_solver
: method to aproximate solutions to the differntail equation. Defaul isTsit5()
.ad_method
:method to evalaute derivatives of the ODE solver. Default isForwardDiffSensitivity()
.
UniversalDiffEq.MultiNODE
— MethodMultiNODE(data;kwargs...)
builds a NODE model to fit to the data. data
is a DataFrame object with time arguments placed in a column labed t
and a second column with a unique index for each time series. The remaining columns have observations of the state variables at each point in time and for each time series.
kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.series_column_name
: Name of column indata
that corresponds to series. Default is"series"
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process erroromega_{proc}
. Default is1.0
.obs_weight
: Weight of observation erroromega_{obs}
. Default is1.0
.reg_weight
: Weight of regularization erroromega_{reg}
. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.ode_solver
: method to aproximate solutions to the differntail equation. Defaul isTsit5()
.ad_method
:method to evalaute derivatives of the ODE solver. Default isForwardDiffSensitivity()
.
UniversalDiffEq.NNDE
— MethodNNDE(data;kwargs ...)
Constructs a nonparametric discrete-time model for the data set data
using a single layer neural network to represent the system's dynamics.
# kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.
UniversalDiffEq.NODE
— MethodNODE(data,X;kwargs ... )
When a data frame X
is supplied the model will run with covariates. The argument X
should have a column for time t
with the value for time in the remaining columns. The values in X
will be interpolated with a linear spline for values of time not included in the data frame.
time_column_name
: Name of column indata
andX
that corresponds to time. Default is"time"
.variable_column_name
: Name of column inX
that corresponds to the variables. Default isnothing
.value_column_name
: Name of column inX
that corresponds to the covariates. Default isnothing
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.ode_solver
: method to aproximate solutions to the differntail equation. Defaul isTsit5()
.ad_method
:method to evalaute derivatives of the ODE solver. Default isForwardDiffSensitivity()
.
UniversalDiffEq.NODE
— MethodNODE(data;kwargs ... )
Constructs a nonparametric continuous-time model for the data set data
using a single layer neural network to represent the system's dynamics.
# kwargs
time_column_name
: Name of column indata
that corresponds to time. Default is"time"
.hidden_units
: Number of neurons in hidden layer. Default is10
.seed
: Fixed random seed for repeatable results. Default is1
.proc_weight
: Weight of process error $omega_{proc}$. Default is1.0
.obs_weight
: Weight of observation error $omega_{obs}$. Default is1.0
.reg_weight
: Weight of regularization error $omega_{reg}$. Default is10^-6
.reg_type
: Type of regularization, whether"L1"
or"L2"
regularization. Default is"L2"
.l
: Extrapolation parameter for forecasting. Default is0.25
.extrap_rho
: Extrapolation parameter for forecasting. Default is0.0
.ode_solver
: method to aproximate solutions to the differntail equation. Defaul isTsit5()
.ad_method
:method to evalaute derivatives of the ODE solver. Default isForwardDiffSensitivity()
.
UniversalDiffEq.NUTS!
— Method NUTS!(UDE, kwargs ...)
Performs Bayesian estimation on the parameters of a UDE using the No U-Turn Sampler (NUTS) algorithm.
kwargs
delta
: Step size used in NUTS adaptor. Default is0.45
.samples
: Number of parameters sampled. Default is500
.burnin
: Number of samples used as burn-in of Bayesian algorithm. Default issamples/10
.verbose
: Should the training loss values be printed? Default istrue
.
UniversalDiffEq.SGLD!
— Method SGLD!(UDE, kwargs ...)
Performs Bayesian estimation on the parameters of an UDE using the Stochastic Gradient Langevin Dynamics (SGLD) sampling algorithm. At each step t
, the stochastic update is provided by a random variable ε with mean 0 and noise η.
math ϵ = a*(b + t-1)^-γ
kwargs
a
: Default is10.0
.b
: Default is1000
.γ
: Default is 0.9.samples
: Number of parameters sampled. Default is500
.burnin
: Number of samples used as burn-in of Bayesian algorithm. Default issamples/10
.verbose
: Should the training loss values be printed? Default istrue
.
UniversalDiffEq.bifurcation_data
— Methodbifurcation_data(model::MultiUDE;N=25)
Calcualtes the equilibrium values of the state variabels $y_t$ as a function of the covariates X_t
and return the value in a data frame. The funciton calcualtes the equilibrium values on a grid of $N$ evenly spaced point for each covariate. The calcualtion are repeated for each time series $i$ included in the training data set.
UniversalDiffEq.bifurcation_data
— Methodbifurcation_data(model::UDE;N=25)
Calcualtes the equilibrium values of the state variabels $y_t$ as a function of the covariates X_t
and return the value in a data frame. The funciton calcualtes the equilibrium values on a grid of $N$ evenly spaced point for each covariate.
UniversalDiffEq.equilibrium_and_stability
— Methodequilibrium_and_stability(UDE,X,lower,upper;t=0,Ntrials=100,tol=10^-3)
Attempts to find all the equilibrium points for the UDE model between the upper and lower bounds, and then returns the real component of the dominant eigenvalue to analyze stability.
...
kwargs
- t = 0: The point in time where the UDE model is evaluated, only relevant for time aware UDEs.
- Ntrials = 100: the number of initializations of the root finding algorithm.
- tol = 10^-3: The threshold Euclidean distance between points beyond which a new equilibrium is sufficiently different to be retained.
...
UniversalDiffEq.equilibrium_and_stability
— Methodequilibrium_and_stability(UDE,lower,upper;t=0,Ntrials=100,tol=10^-3)
Attempts to find all the equilibrium points for the UDE model between the upper and lower bounds, and then returns the real component of the dominant eigenvalue to analyze stability.
...
kwargs
- t = 0: The point in time where the UDE model is evaluated, only relevant for time aware UDEs.
- Ntrials = 100: the number of initializations of the root finding algorithm.
- tol = 10^-3: The threshold Euclidean distance between points beyond which a new equilibrium is sufficiently different to be retained.
...
UniversalDiffEq.equilibrium_and_stability
— Methodequilibrium_and_stability(UDE::MultiUDE,site,X,lower,upper;t=0,Ntrials=100,tol=10^-3)
Attempts to find all the equilibrium points for the UDE model between the upper and lower bounds, and then returns the real component of the dominant eigenvalue to analyze stability.
...
kwargs
- t = 0: The point in time where the UDE model is evaluated, only relevant for time aware UDEs.
- Ntrials = 100: the number of initializations of the root finding algorithm.
- tol = 10^-3: The threshold Euclidean distance between points beyond which a new equilibrium is sufficiently different to be retained.
...
UniversalDiffEq.forecast
— Methodforecast(UDE::BayesianUDE, u0::AbstractVector{}, times::AbstractVector{}; summarize = true, ci = 95)
predictions from the trained model UDE
starting at u0
saving values at times
at each individual sampled parameter. Assumes u0
is the value at time times[1]
If summarize
is true
, this function returns the median prediction as well as the ci
% lower and upper confidence intervals. Othwerise, it returns all the individual predictions for each sampled parameter.
UniversalDiffEq.forecast
— Methodforecast(UDE::BayesianUDE, u0::AbstractVector{}, t0::Real, times::AbstractVector{}; summarize = true, ci = 95)
predictions from the trained model UDE
starting at u0
saving values at times
at each individual sampled parameter. Assumes u0
occurs at time t0
and times
are all larger than t0
.
If summarize
is true
, this function returns the median prediction as well as the ci
% lower and upper confidence intervals. Othwerise, it returns all the individual predictions for each sampled parameter.
UniversalDiffEq.forecast
— Methodforecast(UDE::UDE, u0::AbstractVector{}, times::AbstractVector{})
Predictions from the trained UDE model starting at u0
and saving values at times
. Assumes u0
is the value at initial time times[1]
UniversalDiffEq.forecast
— Methodforecast(UDE::UDE, u0::AbstractVector{}, t0::Real, times::AbstractVector{})
predictions from the trained model UDE
starting at u0
saving values at times
. Assumes u0
occurs at time t0
and times
are all larger than t0
.
UniversalDiffEq.get_NN_parameters
— Methodget_NN_parameters(UDE::UDE)
Returns the values of the weights and biases of the neural network.
UniversalDiffEq.get_parameters
— Methodget_parameters(UDE::UDE)
Returns model parameters.
UniversalDiffEq.get_right_hand_side
— Methodget_right_hand_side(UDE::UDE)
Returns the right-hand side of the differential equation (or difference equation) used to build the process model.
The function will take the state vector u
and time t
if the model does not include covariates. If covariates are included, then the arguments are the state vector u
, covariates vector x
, and time t
.
UniversalDiffEq.gradient_descent!
— Method gradient_descent!(UDE, kwargs ...)
Minimizes the loss function of the UDE
model with the gradient descent algorithm with a step size of step_size
and a maximum number of iterations of maxiter
. Prints the value of the loss function after each iteration when maxiter
is true.
kwargs
step_size
: Step size for ADAM optimizer. Default is0.05
.maxiter
: Maximum number of iterations in gradient descent algorithm. Default is500
.verbose
: Should the training loss values be printed? Default isfalse
.
UniversalDiffEq.leave_future_out
— Methodleave_future_out(model::UDE, training!, k; kwargs... )
Runs a leave future out cross validation on the UDe model model
using the training routine train!
with k
folds. if a path to a csv file is provided usign the path key word then the raw testing data and forecasts will be saved for each fold.
The funtion returns three data frames. The first contains an estimate of the mean aboslue error of the forecasts and assocaited standard error as a fuction of the forecast horizon (1 to k time steps into the future). The second and third are returned in a named tuple with two elements horizon_by_var
and raw
. The data frame horizon_by_var
is containds the forecasting errors seperated by variable and the data frame raw
contains the raw testing and forecasting data. If the model is trained on multiple time series the named tupe will include a third data frame horizon_by_var_by_series
.
# kwargs
path
: the path to a directory to save the output dataframes, defaults to Null
UniversalDiffEq.observation_error_correlations
— Methodobservation_error_correlations(UDE)
The first differnce plot of the observation errors $psilon_t$. This allows the user to check for autocorrelation in the observation errors.
UniversalDiffEq.phase_plane
— Methodphase_plane(UDE::UDE, u0s::AbstractArray; idx=[1,2],T = 100)
Plots the trajectory of state variables as forecasted by the model. Runs a forecast for each provided initial condition out to T
timesteps. Change the state variables that are plotted by changing idx
such that it equals the indexes of the desired state variables as they appear in the data.
UniversalDiffEq.phase_plane
— Methodphase_plane(UDE::UDE; idx=[1,2], u1s=-5:0.25:5, u2s=-5:0.25:5, u3s = 0, T = 100)
Plots the trajectory of state variables as forecasted by the model. Runs a forecast for each permutation of u1 and u2 out to T timesteps. Change the state variables that are plotted by changing idx
such that it equals the indexes of the desired state variables as they appear in the data.
UniversalDiffEq.phase_plane_3d
— Methodphase_plane_3d(UDE::UDE; idx=[1,2,3], u1s=-5:0.25:5, u2s=-5:0.25:5, u3s=-5:0.25:5, T = 100)
The same as phase_plane(), but displays three dimensions/state variables instead of two.
UniversalDiffEq.plot_bifurcation_diagram
— Methodplot_bifurcation_diagram(model::UDE, xvariable; N = 25, color_variable= nothing, conditional_variable = nothing, size= (600, 400))
This function returns a plot of the equilibrium values of the state varaibles $y_t$ as a funciton of the covariates $X_t$. The arguemnt xvariable
determines the covariate plotted on the x-axis. Additional variables can be visualized in sperate panel by specifying the conditional_variable
key word argument or visualized by the color scheme using the color_variable
argument.
The time sereis are treated as an additional covariate that can be visualized by setting the color_variable
or conditional_variable
equal to "series" or the series column name in the training data.
The key word arguent size
controls the dimensions of the final plot.
UniversalDiffEq.plot_bifurcation_diagram
— Methodplot_bifurcation_diagram(model::UDE, xvariable; N = 25, color_variable= nothing, conditional_variable = nothing, size= (600, 400))
This function returns a plot of the equilibrium values of the state varaibles $y_t$ as a funciton of the covariates $X_t$. The arguemnt xvariable
determines the covariate plotted on the x-axis. Additional variables can be visualized in sperate panel by specifying the conditional_variable
key word argument or visualized by the color scheme using the color_variable
argument.
The key word arguent size
controls the dimensions of the final plot.
UniversalDiffEq.plot_forecast
— Methodplot_forecast(UDE::BayesianUDE, test_data::DataFrame)
Plots the model's forecast over the range of the test_data
along with the value of the test data including the median prediction as well as the ci
% lower and upper confidence intervals.
UniversalDiffEq.plot_forecast
— Methodplot_forecast(UDE::BayesianUDE, T::Int)
Plots the models forecast up to T time steps into the future from the last observation including the median prediction as well as the ci
% lower and upper confidence intervals.
UniversalDiffEq.plot_forecast
— Methodplot_forecast(UDE::UDE, test_data::DataFrame)
Plots the model's forecast over the range of the test data along with the value of the test data.
UniversalDiffEq.plot_forecast
— Methodplot_forecast(UDE::UDE, T::Int)
Plots the model's forecast up to T time steps into the future from the last observation.
UniversalDiffEq.plot_predictions
— Methodplot_predictions(UDE::BayesianUDE;ci=95)
Plots the correspondence between the observed state transitons and the predicitons for the model UDE
with a confidence interval ci
.
UniversalDiffEq.plot_predictions
— Methodplot_predictions(UDE::UDE, test_data::DataFrame)
Plots the correspondence between the observed state transitions in test data and the predictions from the UDE model.
UniversalDiffEq.plot_predictions
— Methodplot_predictions(UDE::UDE)
Plots the correspondence between the observed state transitions and the predictions from the UDE model.
UniversalDiffEq.plot_state_estimates
— Methodplot_state_estimates(UDE::UDE)
Plots the values of the state variables estimated by the UDE model.
UniversalDiffEq.predict
— Methodpredict(UDE::BayesianUDE,test_data::DataFrame;summarize = true,ci = 95,df = true)
Uses the Bayesian UDE UDE
to predict the state of the data test_data
for each of the sampled parameters in training.
If summarize
is true
, this function returns the median prediction as well as the ci
% lower and upper confidence intervals. Othwerise, it returns all the individual predictions for each sampled parameter.
If df
is true, this function returns a DataFrame
object. Otherwise, it returns an Array
with the predictions.
UniversalDiffEq.print_parameter_estimates
— Methodprint_parameter_estimates(UDE::UDE)
Prints the values of the known dynamics parameters estimated by the UDE model.
UniversalDiffEq.simulate_coral_data
— Methodsimulate_coral_data()
Generates synthetic coral reef time series data used in the time dependent UDE model example.
UniversalDiffEq.train!
— Methodtrain!(UDE::UDE; kwargs...)
This function provides access to several training routines for UDE models. The user provides the UDE model object and can then choose between several loss functions and optimization algorithms using the keyword arguments. The training routine will update the UDE object with the trained parameters and return any other useful quantities estimated during the training procedure. The default optimizer trains the UDE model by minimizing the loss function using the ADAM gradient descent algorithm for maxiter
steps of size step_size
. Five loss functions are available using the loss_function
argument: conditional likelihood, marginal likelihood, derivative matching, shooting, and multiple shooting. The loss_options
and optim_options
arguments are named tuples that can be used to pass parameters to the training routine.
kwargs
loss_function
: Determines the loss function used to train the model and defaults to "derivative matching" for efficiency.verbose
: If true, the value of the loss function will print between each step of the optimizer.regularization_weight
: Weight given to regularization in the loss function. Default is 0.optimizer
: Determines the optimization algorithm used to train the model and defaults to "ADAM" for gradient descent.loss_options
: Named tuple with keyword arguments to help construct the loss function.optim_options
: Named tuple with keyword arguments to pass to the optimization algorithm. For ADAM, these aremaxiter
, the number of iterations used to run the algorithm, andstep_size
, the size of each iteration.
Loss Functions
Users can choose from one of five loss functions: conditional likelihood, marginal likelihood, derivative matching, shooting, and multiple shooting.
Loss Function | Discrete Model | Continuous Model | Speed |
---|---|---|---|
Conditional likelihood | Yes | Yes | Moderate |
Marginal likelihood | Yes | Yes | Slow |
Derivative matching | No | Yes | Fast |
Shooting | No | Yes | Moderate |
Multiple shooting | No | Yes | Moderate |
Conditional likelihood:
To use the conditional likelihood set the keyword argument loss_function = "conditional likelihood"
.
This option trains the UDE model while accounting for imperfect observations and process uncertainty by maximizing the conditional likelihood of a state-space model where the UDE is used as the process model. The conditional likelihood is faster to compute, but can be less accurate than the marginal likelihood.
Marginal likelihood:
To use the marginal likelihood set the keyword argument loss_function = "marginal likelihood"
.
This option maximizes the marginal likelihood of a state-space model, which is approximated using an unscented Kalman filter. This option is slower than the conditional likelihood but should, in theory, increase the accuracy of the trained model (i.e., reduce bias).
loss_options
process_error
: An initial estimate of the level of process error. Default is 0.1.observation_error
: The level of observation error in the data set. There is no default, so it will throw an error if not provided.α
: Parameter for the Kalman filter algorithm. Default is 10^-3.β
: Parameter for the Kalman filter algorithm. Default is 2.κ
: Parameter for the Kalman filter algorithm. Default is 0.
Derivative matching:
To use the derivative matching training routine set loss_function = "derivative matching"
.
This function trains the UDE model in a two-step process. First, a smoothing function is fit to the data using a spline regression. Then, the UDE model is trained by comparing the derivatives of the smoothing functions to the derivatives predicted by the right-hand side of the UDE. This training routine is much faster than the alternatives, but may be less accurate.
loss_options
d
: The number of degrees of freedom in the curve fitting function. Defaults to 12.alg
: The algorithm used to fit the curve to the data set. See the DataInterpolations package for details. Defaults to generalized cross-validation:gcv_svd
.remove_ends
: The number of data points to leave off of the end of the data set when training the UDE to reduce edge effects from the curve fitting process. Defaults to 0.
Shooting:
To use the shooting training routine set loss_function = "shooting"
.
This option calculates the loss by solving the ODE from the initial to the final data point and comparing the observed to the predicted trajectory with mean squared error (MSE). The initial data point is estimated as a free parameter to reduce the impacts of observation error.
Multiple shooting:
To use the multiple shooting training routine set loss_function = "multiple shooting"
.
This option calculates the loss by breaking the data into blocks of sequential observations. It then uses the UDE to forecast from the initial data point in each block to the first data point in the next block. The loss is defined as the MSE between the forecasts and the data points. The initial data point in each block is estimated as a free parameter to reduce the impacts of observation error.
loss_options
pred_length
: The number of data points in each block. The default is 10.
Optimization Algorithms
The method used to minimize the loss function. Two options are available: ADAM and BFGS. ADAM is a first-order gradient descent algorithm, while BFGS is a quasi-Newton method that uses approximate second-order information.
The user can specify the maximum number of iterations maxiter
for each algorithm using the optim_options
keyword argument. For ADAM the optim_options
can be used to specify the step size step_size
and for BFGS you can specify the initial step norm initial_step_norm
.
UniversalDiffEq.vectorfield_and_nullclines
— Methodvectorfieldandnullclines(UDE; kwargs)
Calculate the vector field and nullclines of the 2D UDE
model and returns their plot.
kwargs
-t
: Time step t
at which the vector field and nullclines are calculated. Default is 0
. -n
: Number of elements per axes to evaluate vector field at. Default is 15
. -lower
: Lower limits of vector field and nullclines. Default is [0.0,0.0]
. -upper
: Upper limits of vector field and nullclines. Default is [1.0,1.0]
. -arrowlength
: Arrow size of vector field plot. Default is 2
. -arrow_color
: Arrow color of vector field plot. Default is grey
. -xlabel
: X-label of vector field plot. Default is u1
. -ylabel
: Y-label of vector field plot. Default is u2
. -title
: Plot title. Default is Vector field
. -color_u1
: Color of nullcline in x-axis. Default is "red"
. -color_u2
: Color of nullcline in y-axis. Default is "black"
. -legend
: Position of legends of nullcines in plot. Default is :outerright
.