API Reference {#Reference}
Types
Structs
Main.Sarimax.SARIMAModel
— TypeThe SARIMAModel
struct represents a SARIMA model. It contains the following fields:
y
: The time series data.p
: The autoregressive order for the non-seasonal part.d
: The degree of differencing.q
: The moving average order for the non-seasonal part.seasonality
: The seasonality period.P
: The autoregressive order for the seasonal part.D
: The degree of seasonal differencing.Q
: The moving average order for the seasonal part.metadata
: A dictionary containing model metadata.exog
: Optional exogenous variables.c
: The constant term.trend
: The trend term.ϕ
: The autoregressive coefficients for the non-seasonal part.θ
: The moving average coefficients for the non-seasonal part.Φ
: The autoregressive coefficients for the seasonal part.Θ
: The moving average coefficients for the seasonal part.ϵ
: The residuals.exogCoefficients
: The coefficients of the exogenous variables.σ²
: The variance of the residuals.fitInSample
: The in-sample fit.forecast
: The forecast.silent
: Whether to suppress output.allowMean
: Whether to include a mean term in the model.allowDrift
: Whether to include a drift term in the model.keepProvidedCoefficients
: Whether to keep the provided coefficients.
Enums
Main.Sarimax.Datasets
— TypeThe Datasets Enum is used to identify the dataset used in the loadDataset
function.
The Datasets Enum is defined as follows:
AIR_PASSENGERS = 1
GDPC1 = 2
NROU = 3
The loadDataset
function uses this Enum to determine the dataset to be loaded.
Exceptions
Main.Sarimax.ModelNotFitted
— Typestruct ModelNotFitted <: Exception
An exception type that indicates the model has not been fitted yet.
Usage
This exception can be thrown when an operation that requires a fitted model is attempted on an unfitted model.
Main.Sarimax.MissingMethodImplementation
— TypeMissingMethodImplementation <: Exception
Custom exception type for indicating that a required method is not implemented in the model.
Fields
method::String
: The name of the method that is missing.
Main.Sarimax.MissingExogenousData
— TypeMissingExogenousData <: Exception
An exception type that indicates the absence of exogenous data required for forecasting the requested horizon.
Main.Sarimax.InconsistentDatePattern
— Typestruct InconsistentDatePattern <: Exception
An exception type to indicate that the timestamps do not follow a consistent pattern.
Main.Sarimax.InvalidParametersCombination
— TypeInvalidParametersCombination <: Exception
A custom exception type to indicate that the combination of parameters provided to the model is invalid.
Fields
msg::String
: A message describing why the parameters are invalid.
Constructors
Main.Sarimax.SARIMA
— MethodSARIMA constructor.
Parameters:
- y: TimeArray with the time series.
- p: Int with the autoregressive order for the non-seasonal part.
- d: Int with the degree of differencing.
- q: Int with the moving average order for the non-seasonal part.
- seasonality: Int with the seasonality period.
- P: Int with the autoregressive order for the seasonal part.
- D: Int with the degree of seasonal differencing.
- Q: Int with the moving average order for the seasonal part.
- silent: Bool to supress output.
- allowMean: Bool to include a mean term in the model.
- allowDrift: Bool to include a drift term in the model.
Main.Sarimax.SARIMA
— MethodSARIMA constructor to initialize model with provided coefficients.
Parameters:
- y: TimeArray with the time series.
- exog: TimeArray with the exogenous variables.
- arCoefficients: Vector with the autoregressive coefficients.
- maCoefficients: Vector with the moving average coefficients.
- seasonalARCoefficients: Vector with the autoregressive coefficients for the seasonal component.
- seasonalMACoefficients: Vector with the moving average coefficients for the seasonal component.
- mean: Float with the mean term.
- trend: Float with the trend term.
- exogCoefficients: Vector with the exogenous coefficients.
- d: Int with the degree of differencing.
- D: Int with the degree of seasonal differencing.
- seasonality: Int with the seasonality period.
- silent: Bool to supress output.
- allowMean: Bool to include a mean term in the model.
- allowDrift: Bool to include a drift term in the model.
Main.Sarimax.SARIMA
— MethodSARIMA constructor.
Parameters:
- y: TimeArray with the time series.
- exog: TimeArray with the exogenous variables.
- p: Int with the order of the AR component.
- d: Int with the degree of differencing.
- q: Int with the order of the MA component.
- seasonality: Int with the seasonality period.
- P: Int with the order of the seasonal AR component.
- D: Int with the degree of seasonal differencing.
- Q: Int with the order of the seasonal MA component.
- silent: Bool to supress output.
- allowMean: Bool to include a mean term in the model.
- allowDrift: Bool to include a drift term in the model.
Model Functions
Model Fitting and Prediction
Main.Sarimax.fit!
— Functionfit!(
model::SARIMAModel;
silent::Bool=true,
optimizer::DataType=Ipopt.Optimizer,
objectiveFunction::String="mse"
automaticExogDifferentiation::Bool=false
)
Estimate the SARIMA model parameters via non-linear least squares. The resulting optimal parameters as well as the residuals and the model's σ² are stored within the model. The default objective function used to estimate the parameters is the mean squared error (MSE), but it can be changed to the maximum likelihood (ML) by setting the objectiveFunction
parameter to "ml".
Arguments
model::SARIMAModel
: The SARIMA model to be fitted.silent::Bool
: Whether to suppress solver output. Default istrue
.optimizer::DataType
: The optimizer to be used for optimization. Default isIpopt.Optimizer
.objectiveFunction::String
: The objective function used for estimation. Default is "mse".automaticExogDifferentiation::Bool
: Whether to automatically differentiate the exogenous variables. Default isfalse
.
Example
julia> airPassengers = loadDataset(AIR_PASSENGERS)
julia> model = SARIMA(airPassengers,0,1,1;seasonality=12,P=0,D=1,Q=1)
julia> fit!(model)
Main.Sarimax.predict!
— Functionpredict!(
model::SARIMAModel;
stepsAhead::Int = 1
seed::Int = 1234,
isSimulation::Bool = false,
displayConfidenceIntervals::Bool = false,
confidenceLevel::Fl = 0.95
automaticExogDifferentiation::Bool=false
) where Fl<:AbstractFloat
Predicts the SARIMA model for the next stepsAhead
periods. The resulting forecast is stored within the model in the forecast
field.
Arguments
model::SARIMAModel
: The SARIMA model to make predictions.stepsAhead::Int
: The number of periods ahead to forecast (default: 1).seed::Int
: Seed for random number generation when simulating forecasts (default: 1234).isSimulation::Bool
: Whether to perform a simulation-based forecast (default: false).displayConfidenceIntervals::Bool
: Whether to display confidence intervals (default: false).confidenceLevel::Fl
: The confidence level for the confidence intervals (default: 0.95).automaticExogDifferentiation::Bool
: Whether to automatically differentiate the exogenous variables. Default isfalse
.
Example
```julia julia> airPassengers = loadDataset(AIR_PASSENGERS)
julia> model = SARIMA(airPassengers, 0, 1, 1; seasonality=12, P=0, D=1, Q=1)
julia> fit!(model)
julia> predict!(model; stepsAhead=12)
Main.Sarimax.auto
— Functionauto(
y::TimeArray;
exog::Union{TimeArray,Nothing}=nothing,
seasonality::Int=1,
d::Int = -1,
D::Int = -1,
maxp::Int = 5,
maxd::Int = 2,
maxq::Int = 5,
maxP::Int = 2,
maxD::Int = 1,
maxQ::Int = 2,
maxOrder::Int = 5,
informationCriteria::String = "aicc",
allowMean:Union{Bool,Nothing} = nothing,
allowDrift::Union{Bool,Nothing} = nothing,
integrationTest::String = "kpss",
seasonalIntegrationTest::String = "seas",
objectiveFunction::String = "mse",
assertStationarity::Bool = true,
assertInvertibility::Bool = true,
showLogs::Bool = false,
outlierDetection::Bool = false
searchMethod::String = "stepwise"
)
Automatically fits the best SARIMA model according to the specified parameters.
Arguments
y::TimeArray
: The time series data.exog::Union{TimeArray,Nothing}
: Optional exogenous variables. IfNothing
, no exogenous variables are used.seasonality::Int
: The seasonality period. Default is 1 (non-seasonal).d::Int
: The degree of differencing for the non-seasonal part. Default is -1 (auto-select).D::Int
: The degree of differencing for the seasonal part. Default is -1 (auto-select).maxp::Int
: The maximum autoregressive order for the non-seasonal part. Default is 5.maxd::Int
: The maximum integration order for the non-seasonal part. Default is 2.maxq::Int
: The maximum moving average order for the non-seasonal part. Default is 5.maxP::Int
: The maximum autoregressive order for the seasonal part. Default is 2.maxD::Int
: The maximum integration order for the seasonal part. Default is 1.maxQ::Int
: The maximum moving average order for the seasonal part. Default is 2.maxOrder::Int
: The maximum order for the non-seasonal part. Default is 5.informationCriteria::String
: The information criteria to be used for model selection. Options are "aic", "aicc", or "bic". Default is "aicc".allowMean::Union{Bool,Nothing}
: Whether to include a mean term in the model. Default is nothing.allowDrift::Union{Bool,Nothing}
: Whether to include a drift term in the model. Default is nothing.integrationTest::String
: The integration test to be used for determining the non-seasonal integration order. Default is "kpss".seasonalIntegrationTest::String
: The integration test to be used for determining the seasonal integration order. Default is "seas".objectiveFunction::String
: The objective function to be used for model selection. Options are "mse", "ml", or "bilevel". Default is "mse".assertStationarity::Bool
: Whether to assert stationarity of the fitted model. Default is true.assertInvertibility::Bool
: Whether to assert invertibility of the fitted model. Default is true.showLogs::Bool
: Whether to suppress output. Default is false.outlierDetection::Bool
: Whether to perform outlier detection. Default is false.- searchMethod::String = "stepwise"
References
- Hyndman, RJ and Khandakar. "Automatic time series forecasting: The forecast package for R." Journal of Statistical Software, 26(3), 2008.
Main.Sarimax.simulate
— Functionsimulate(
model::SARIMAModel,
stepsAhead::Int = 1,
numScenarios::Int = 200,
seed::Int = 1234
)
Simulates the SARIMA model for the next stepsAhead
periods assuming that the model's estimated σ². Returns a vector of numScenarios
scenarios of the forecasted values.
Arguments
model::SARIMAModel
: The SARIMA model to simulate.stepsAhead::Int
: The number of periods ahead to simulate. Default is 1.numScenarios::Int
: The number of simulation scenarios. Default is 200.seed::Int
: The seed of the simulation. Default is 1234.
Returns
Vector{Vector{AbstractFloat}}
: A vector of scenarios, each containing the forecasted values for the nextstepsAhead
periods.
Example
julia> airPassengers = loadDataset(AIR_PASSENGERS)
julia> model = SARIMA(airPassengers, 0, 1, 1; seasonality=12, P=0, D=1, Q=1)
julia> fit!(model)
julia> scenarios = simulate(model, stepsAhead=12, numScenarios=1000)
Model Evaluation
Main.Sarimax.loglikelihood
— Functionloglikelihood(model::SarimaxModel)
Calculate the log-likelihood of a SARIMAModel. The log-likelihood is calculated based on the formula -0.5 * (T * log(2π) + T * log(σ²) + sum(ϵ.^2 ./ σ²))
where:
- T is the length of the residuals vector (ϵ).
- σ² is the estimated variance of the model.
Arguments
model::SarimaxModel
: A SARIMAModel object.
Returns
- The log-likelihood of the SARIMAModel.
Errors
MissingMethodImplementation("fit!")
: Thrown if thefit!
method is not implemented for the given model type.ModelNotFitted()
: Thrown if the model has not been fitted.
Main.Sarimax.loglike
— Functionloglike(model::SarimaxModel)
Calculate the log-likelihood of a SARIMAModel using the normal probability density function. The log-likelihood is calculated by summing the logarithm of the probability density function (PDF) of each data point under the assumption of a normal distribution with mean 0 and standard deviation equal to the square root of the estimated variance (σ²) of the model.
Arguments
model::SarimaxModel
: A SARIMAModel object.
Returns
- The log-likelihood of the SARIMAModel.
Errors
MissingMethodImplementation("fit!")
: Thrown if thefit!
method is not implemented for the given model type.ModelNotFitted()
: Thrown if the model has not been fitted.
Main.Sarimax.aic
— Functionaic(K::Int, loglikeVal::Fl) where Fl<:AbstractFloat -> Fl
Calculate the Akaike Information Criterion (AIC) for a given number of parameters and log-likelihood value.
Arguments
K::Int
: Number of parameters in the model.loglikeVal::Fl
: Log-likelihood value of the model.
Returns
The AIC value calculated using the formula: AIC = 2K - 2loglikeVal.
aic(model::SarimaxModel; offset::Fl) -> Fl where Fl<:AbstractFloat
Calculate the Akaike Information Criterion (AIC) for a Sarimax model.
Arguments
model::SarimaxModel
: The Sarimax model for which AIC is calculated.offset::Fl=0.0
: Offset value to be added to the AIC value.
Returns
The AIC value calculated using the number of parameters and log-likelihood value of the model.
Errors
- Throws a
MissingMethodImplementation
if thegetHyperparametersNumber
method is not implemented for the given model type.
Main.Sarimax.aicc
— Functionaicc(T::Int, K::Int, loglikeVal::Fl) where Fl<:AbstractFloat -> Fl
Calculate the corrected Akaike Information Criterion (AICc) for a given number of observations, number of parameters, and log-likelihood value.
Arguments
T::Int
: Number of observations in the data.K::Int
: Number of parameters in the model.loglikeVal::Fl
: Log-likelihood value of the model.
Returns
The AICc value calculated using the formula: AICc = AIC(K, loglikeVal) + ((2KK + 2*K) / (T - K - 1)).
aicc(model::SarimaxModel; offset::Fl) -> Fl where Fl<:AbstractFloat
Calculate the Corrected Akaike Information Criterion (AICc) for a Sarimax model.
Arguments
model::SarimaxModel
: The Sarimax model for which AICc is calculated.offset::Fl=0.0
: Offset value to be added to the AICc value.
Returns
The AICc value calculated using the number of parameters, sample size, and log-likelihood value of the model.
Errors
- Throws a
MissingMethodImplementation
if thegetHyperparametersNumber
method is not implemented for the given model type.
Main.Sarimax.bic
— Functionbic(T::Int, K::Int, loglikeVal::Fl) -> Fl
Calculate the Bayesian Information Criterion (BIC) for a given number of observations, number of parameters, and log-likelihood value.
Arguments
T::Int
: Number of observations in the data.K::Int
: Number of parameters in the model.loglikeVal::Fl
: Log-likelihood value of the model.
Returns
The BIC value calculated using the formula: BIC = log(T) * K - 2 * loglikeVal.
bic(model::SarimaxModel; offset::Fl) -> Fl where Fl<:AbstractFloat
Calculate the Bayesian Information Criterion (BIC) for a Sarimax model.
Arguments
model::SarimaxModel
: The Sarimax model for which BIC is calculated.offset::Fl=0.0
: Offset value to be added to the BIC value.
Returns
The BIC value calculated using the number of parameters, sample size, and log-likelihood value of the model.
Errors
- Throws a
MissingMethodImplementation
if thegetHyperparametersNumber
method is not implemented for the given model type.
Time Series Operations
Main.Sarimax.differentiate
— Functiondifferentiate(
series::TimeArray,
d::Int=0,
D::Int=0,
s::Int=1
)
Differentiates a TimeArray
series
d
times and D
times with a seasonal difference of s
periods.
Arguments
series::TimeArray
: The time series data to differentiate.d::Int=0
: The number of non-seasonal differences to take.D::Int=0
: The number of seasonal differences to take.s::Int=1
: The seasonal period for the differences.
Returns
A differentiated TimeArray
.
Errors
- This method only works with
d
andD
in the set {0,1}.
Example
julia> airPassengers = loadDataset(AIR_PASSENGERS)
julia> stationaryAirPassengers = differentiate(airPassengers, d=1, D=1, s=12)
Main.Sarimax.integrate
— Functionintegrate(initialValues::Vector{Fl}, diffSeries::Vector{Fl}, d::Int, D::Int, s::Int) where Fl<:AbstractFloat
Converts a differentiated time series back to its original scale.
Arguments
initialValues::Vector{Fl}
: Initial values of the original time series.diffSeries::Vector{Fl}
: Differentiated time series.d::Int
: Order of non-seasonal differencing.D::Int
: Order of seasonal differencing.s::Int
: Seasonal period.
Returns
origSeries::Vector{Fl}
: Time series in the original scale.
Main.Sarimax.differentiatedCoefficients
— FunctiondifferentiatedCoefficients(d::Int, D::Int, s::Int, Fl::DataType=Float64)
Compute the coefficients for differentiating a time series.
Arguments
d::Int
: Order of non-seasonal differencing.D::Int
: Order of seasonal differencing.s::Int
: Seasonal period.Fl
: The type of the coefficients. Default isFloat64
.
Returns
coeffs::Vector{AbstractFloat}
: Coefficients for differentiation.
Main.Sarimax.toMA
— FunctiontoMA(model::SARIMAModel, maxLags::Int=12)
Convert a SARIMA model to a Moving Average (MA) model.
# Arguments
- `model::SARIMAModel`: The SARIMA model to convert.
- `maxLags::Int=12`: The maximum number of lags to include in the MA model.
# Returns
- `MAmodel::MAModel`: The coefficients of the lagged errors in the MA model.
# References
- Brockwell, P. J., & Davis, R. A. Time Series: Theory and Methods (page 92). Springer(2009)
Dataset and Utility Functions
Main.Sarimax.loadDataset
— FunctionloadDataset(
dataset::Datasets
)
Loads a dataset from the Datasets
enum.
Example
julia> airPassengers = loadDataset(AIR_PASSENGERS)
204×1 TimeArray{Float64, 1, Date, Vector{Float64}} 1991-07-01 to 2008-06-01
│ │ value │
├────────────┼─────────┤
│ 1991-07-01 │ 3.5266 │
│ 1991-08-01 │ 3.1809 │
│ ⋮ │ ⋮ │
│ 2008-06-01 │ 19.4317 │
loadDataset(
df::DataFrame,
showLogs::Bool=false
)
Loads a dataset from a Dataframe. If the DataFrame does not have a column named date
a new column will be created with the index of the DataFrame.
Arguments
df::DataFrame
: The DataFrame to be converted to a TimeArray.showLogs::Bool=false
: If true, logs will be shown.
Example
julia> airPassengersDf = CSV.File("datasets/airpassengers.csv") |> DataFrame
julia> airPassengers = loadDataset(airPassengersDf)
204×1 TimeArray{Float64, 1, Date, Vector{Float64}} 1991-07-01 to 2008-06-01
│ │ value │
├────────────┼─────────┤
│ 1991-07-01 │ 3.5266 │
│ 1991-08-01 │ 3.1809 │
│ ⋮ │ ⋮ │
│ 2008-06-01 │ 19.4317 │
Main.Sarimax.splitTrainTest
— FunctionsplitTrainTest(
data::TimeArray;
trainSize::Fl=0.8
) where Fl<:AbstractFloat
Splits the time series in training and testing sets.
Main.Sarimax.identifyGranularity
— FunctionidentifyGranularity(datetimes::Vector{T})
Identifies the granularity of an array of timestamps.
Arguments
datetimes::Vector{T}
: An array of TimeType objects.
Returns
A tuple (granularity, frequency, weekdays)
where:
granularity
: The identified granularity, which could be one of [:Second, :Minute, :Hour, :Day, :Week, :Month, :Year].frequency
: The frequency of the identified granularity.weekdays
: A boolean indicating whether the series uses weekdays only.
Errors
Throws an error if the timestamps do not follow a consistent pattern.
Main.Sarimax.buildDatetimes
— FunctionbuildDatetimes(startDatetime, granularity, weekDaysOnly, datetimesLength)
Builds an array of DateTime objects based on a given starting DateTime, granularity, and length.
Arguments
startDatetime::T
: The DateTime from which the dateTime array will be computed. It won't be included in the final arraygranularity::Dates.Period
: The granularity by which to increment the timestamps.weekDaysOnly::Bool
: Whether to include only weekdays (Monday to Friday) in the timestamps.datetimesLength::Int
: The length of the array of DateTime objects to build.
Returns
An array of DateTime objects.
Model Information
Main.Sarimax.hasFitMethods
— FunctionhasFitMethods(modelType::Type{<:SarimaxModel}) -> Bool
Check if a given SarimaxModel
type has the fit!
method implemented.
Arguments
modelType::Type{<:SarimaxModel}
: Type of the Sarimax model to check.
Returns
A boolean indicating whether the fit!
method is implemented for the specified model type.
Main.Sarimax.hasHyperparametersMethods
— FunctionhasHyperparametersMethods(modelType::Type{<:SarimaxModel}) -> Bool
Checks if a given SarimaxModel
type has methods related to hyperparameters.
Arguments
modelType::Type{<:SarimaxModel}
: Type of the Sarimax model to check.
Returns
A boolean indicating whether the hyperparameter-related methods are implemented for the specified model type.
Main.Sarimax.getHyperparametersNumber
— FunctiongetHyperparametersNumber(model::SARIMAModel)
Returns the number of hyperparameters of a SARIMA model.
Arguments
model::SARIMAModel
: The SARIMA model.
Returns
Int
: The number of hyperparameters.
Model Manipulation
Base.print
— Functionprint(model::SARIMAModel)
Prints the SARIMA model.
Arguments
model::SARIMAModel
: The SARIMA model to print.
Example
julia> model = SARIMA(1, 0, 1; P=1, D=0, Q=1, seasonality=12, allowMean=true, allowDrift=false)
julia> print(model)
Main.Sarimax.copyTimeArray
— MethodcopyTimeArray(y::TimeSeries.TimeArray)
Create a shallow copy of a TimeArray object.
This function creates a new TimeArray with copies of the timestamp and values from the original TimeArray. The new TimeArray is independent of the original, but the elements themselves are not deeply copied.
Arguments
y::TimeSeries.TimeArray
: The TimeArray to copy.
Returns
TimeSeries.TimeArray
: A new TimeArray with copies of the timestamp and values.
Examples
original = TimeArray(Date.(2021:2023), [1, 2, 3])
copied = copyTimeArray(original)
Main.Sarimax.deepcopyTimeArray
— MethoddeepcopyTimeArray(y::TimeSeries.TimeArray)
Create a deep copy of a TimeArray object.
This function creates a new TimeArray with deep copies of the timestamp and values from the original TimeArray. The new TimeArray and all its elements are completely independent of the original.
Arguments
y::TimeSeries.TimeArray
: The TimeArray to deep copy.
Returns
TimeSeries.TimeArray
: A new TimeArray with deep copies of the timestamp and values.
Examples
original = TimeArray(Date.(2021:2023), [1, 2, 3])
deepCopied = deepcopyTimeArray(original)