API Reference

This page provides documentation for all exported functions and types in NowcastAutoGP.jl.

Index

API

NowcastAutoGP.IntegratedBrownianMotionType
IntegratedBrownianMotion(origin[, amplitude=1])

Once-integrated Brownian motion (integrated Wiener process) covariance kernel.

\[k(t, t') = \theta_2 \, \frac{a^2 (3b - a)}{6}, \qquad a = \min(t, t') - \theta_1, \quad b = \max(t, t') - \theta_1\]

This is the covariance of $X(t) = \int_{\theta_1}^{t} W(s)\,\mathrm{d}s$, the integral of a Brownian motion that starts at the origin $\theta_1$ (where both the value and the variance are zero); amplitude $\theta_2$ scales the variance. Draws are smoother (once-differentiable) than RandomWalk draws — an integrated-random-walk prior, a natural choice for trends whose rate of change drifts like a random walk.

The kernel is positive-definite only when origin <= min(t) over the evaluated time points, so that every a, b >= 0.

Like RandomWalk, this kernel is defined in NowcastAutoGP (not AutoGP) by extending AutoGP's GP interface (eval_cov, reparameterize, rescale).

Connection to cubic splines

The predictive mean of a GP with kernel composition of primitives Const + Linear + IntegratedBrownianMotion is a cubic spline (piecewise cubic polynomial) between data points and linear outside the data range, see [Rasmussen & Williams, 2006, §6.3]. This means that a GP with this kernel composition represents a Bayesian cubic spline model, with the IntegratedBrownianMotion kernel controlling the smoothness of the spline. Note that the kernel covariance here seems to differ in form from the one in Rasmussen & Williams, but noting that |t - t'| = b - a, the two forms are equivalent.

source
NowcastAutoGP.RandomWalkType
RandomWalk(origin[, amplitude=1])

Random walk (Wiener process) covariance kernel.

\[k(t, t') = \theta_2 \left( \min(t, t') - \theta_1 \right)\]

The origin $\theta_1$ is the time at which the process starts (where the variance is zero); amplitude $\theta_2$ scales the per-unit-time variance. Draws from this kernel are continuous-time random walks (Brownian motion): non-stationary, with variance growing linearly away from origin.

The kernel is positive-definite only when origin <= min(t) over the evaluated time points, so that every min(t, t') - origin >= 0. This is the GP analogue of an i.i.d.-increment random walk.

This kernel is defined in NowcastAutoGP (not AutoGP) by extending AutoGP's GP interface (eval_cov, reparameterize, rescale). It mirrors AutoGP's primitive kernel structure.

source
NowcastAutoGP.TDataType
TData{D, F}

A container for transformed time series data used in nowcasting models.

Type Parameters

  • D: Type for dates/timestamps (e.g., Date, DateTime)
  • F: Type for numeric values, automatically promoted from input types

Fields

  • ds::Vector{D}: Vector of dates or timestamps corresponding to observations
  • y::Vector{F}: Vector of transformed target values (result of applying transformation)
  • values::Vector{F}: Vector of original values, converted to common type with y

Constructor

TData(ds, values; transformation)

Create a TData instance by applying a transformation to the input values.

Arguments

  • ds: Vector of dates or timestamps
  • values: Vector of original numeric values
  • transformation: Function to apply element-wise to values to create y

The constructor automatically promotes types using promote_type to ensure y and values have compatible numeric types.

Example

using Dates

dates = [Date(2023, 1, 1), Date(2023, 1, 2), Date(2023, 1, 3)]
raw_values = [10, 20, 30]

# Apply log transformation
tdata = TData(dates, raw_values; transformation = log)

# Apply custom transformation
tdata = TData(dates, raw_values; transformation = x -> (x - mean(raw_values)) / std(raw_values))

Validation

The constructor ensures that ds and values have the same length and throws an ArgumentError if they don't match.

source
NowcastAutoGP._get_offsetMethod
_get_offset(values::Vector{F}) where {F <: Real}

Internal function to compute an offset for transformations to ensure numerical stability.

source
NowcastAutoGP._inv_boxcoxMethod
_inv_boxcox(λ::Real, offset::F, max_values) where {F}

Internal function to compute the inverse Box-Cox transformation with edge case handling.

source
NowcastAutoGP._stabilize_for_fitMethod
_stabilize_for_fit(y; flat_threshold=1.0e-3)

Guard against a degenerate (near-constant) transformed series before fitting.

Gaussian process fitting needs the data to have some spread: AutoGP rescales y by its range, so a genuinely flat series makes the covariance singular (PosDefException, issue #51). When the relative range of y ((max - min) / (|mean| + 1)) falls below flat_threshold, add small Gaussian jitter so the series is fittable; otherwise return y unchanged. Operates on the (already transformed) values, so it is independent of which transformation was used.

This is a secondary safety net for data that is flat even after transformation. The common case of a pathological Box-Cox λ on otherwise-fittable data is handled earlier, in get_transformations, by falling back to a log transformation.

source
NowcastAutoGP.create_nowcast_dataMethod
create_nowcast_data(nowcasts::AbstractMatrix, dates::Vector{Date}; transformation = y -> y)

Create nowcast data structures from a matrix of nowcast scenarios.

Arguments

  • nowcasts: A matrix where each column represents one nowcast scenario. The number of rows must match the length of dates.
  • dates: A vector of Date objects corresponding to the nowcast time points.
  • transformation: A function to apply to the nowcast values (default: identity).

Returns

A vector of NamedTuples, where each NamedTuple represents one nowcast scenario with fields:

  • ds: The dates vector
  • y: The transformed nowcast values
  • values: The original (untransformed) nowcast values

Notes

This method converts the matrix to a vector of columns internally and delegates to the vector method.

Example

# Matrix with 3 time points (rows) and 2 scenarios (columns)
nowcasts = [10.5 9.8; 11.2 10.9; 12.1 11.5]
dates = [Date(2024,1,1), Date(2024,1,2), Date(2024,1,3)]
nowcast_data = create_nowcast_data(nowcasts, dates; transformation = log)
# Returns vector of 2 NamedTuples, each with transformed and original values
source
NowcastAutoGP.create_nowcast_dataMethod
create_nowcast_data(nowcasts::AbstractVector, dates::Vector{Date}; transformation = y -> y)

Create nowcast data structures from a vector of nowcast scenarios.

Arguments

  • nowcasts: A vector where each element is a vector of nowcast values representing one scenario. All inner vectors must have the same length as dates.
  • dates: A vector of Date objects corresponding to the nowcast time points.
  • transformation: A function to apply to the nowcast values (default: identity).

Returns

A vector of NamedTuples, where each NamedTuple represents one nowcast scenario with fields:

  • ds: The dates vector
  • y: The transformed nowcast values
  • values: The original (untransformed) nowcast values

Example

# Two nowcast scenarios for 3 dates
nowcasts = [[10.5, 11.2, 12.1], [9.8, 10.9, 11.5]]
dates = [Date(2024,1,1), Date(2024,1,2), Date(2024,1,3)]
nowcast_data = create_nowcast_data(nowcasts, dates; transformation = log)
# Returns vector of 2 NamedTuples, each with transformed and original values
source
NowcastAutoGP.forecastMethod
forecast(model, forecast_dates, forecast_draws::Int)

Generate forecast samples from a fitted AutoGP model.

Arguments

  • model: Fitted AutoGP.GPModel.
  • forecast_dates: Vector or range of dates to predict.
  • forecast_draws: Number of samples to draw.

Keyword arguments

  • inv_transformation: Function applied elementwise to map forecasts back to the original scale (default: identity).
  • forecast_n_hmc: If nothing, draw from the current model state. If an Int, run that many HMC parameter steps before each draw (default: nothing).

Returns

  • A matrix of samples with size (length(forecast_dates), forecast_draws).
source
NowcastAutoGP.forecast_with_nowcastsMethod
forecast_with_nowcasts(base_model, nowcasts, forecast_dates, forecast_draws_per_nowcast;
                      inv_transformation = y -> y, n_mcmc = 0, n_hmc = 0, ess_threshold = 0.0)

Generate forecasts by conditioning on multiple nowcast scenarios.

Arguments

  • base_model: Fitted AutoGP.GPModel trained on confirmed (non-nowcast) data.
  • nowcasts: Vector of TData scenarios with fields ds, y, and values.
  • forecast_dates: Vector or range of dates to predict.
  • forecast_draws_per_nowcast: Samples per scenario.

Keyword arguments

  • inv_transformation: Function applied elementwise to map forecasts back to the original scale (default: identity).
  • n_mcmc: Number of MCMC structure steps after adding each nowcast (default: 0). If > 0, n_hmc must also be > 0.
  • n_hmc: Number of HMC parameter steps per MCMC step (default: 0). Can be > 0 even if n_mcmc == 0.
  • ess_threshold: Effective sample size threshold for particle resampling, as a fraction of total particles (default: 0.0).
  • forecast_n_hmc: Number of HMC steps to run before each forecast draw (default: nothing). If nothing, no HMC steps are taken during forecasting.
  • verbose: If true, display progress information during forecasting (default: false).

Returns

  • A matrix with size (length(forecast_dates), length(nowcasts) * forecast_draws_per_nowcast).

Notes

  • Each scenario operates on an independent copy of the base model, so the original model is never mutated.
  • n_mcmc == 0 && n_hmc > 0 performs parameter-only updates to the particle ensemble; n_mcmc > 0 && n_hmc > 0 performs full MCMC.
  • forecast_n_hmc is independent of n_mcmc and n_hmc and controls HMC steps only during forecasting, not during nowcast incorporation.

If n_mcmc == 0 && n_hmc == 0 && forecast_n_hmc > 0, HMC steps are only taken during forecasting, not during nowcast incorporation.

Example

nowcast_scenarios = [
    (ds = [Date(2024,1,1), Date(2024,1,2)], y = [10.5, 11.2], values = [10.5, 11.2]),
    (ds = [Date(2024,1,1), Date(2024,1,2)], y = [9.8, 10.9], values = [9.8, 10.9]),
]
forecast_dates = Date(2024,1,1):Day(1):Date(2024,1,10)
forecasts = forecast_with_nowcasts(base_model, nowcast_scenarios, forecast_dates, 100)
source
NowcastAutoGP.get_transformationsMethod
get_transformations(transform_name::String, values::Vector{F}) where {F <: Real}

Return a tuple of transformation and inverse transformation functions for the specified transformation type.

This function creates appropriate data transformations for Gaussian Process modeling, where the goal is to transform the input data to make it more suitable for modeling (typically more Gaussian-like) and then provide the inverse transformation to convert predictions back to the original scale.

Arguments

  • transform_name::String: The name of the transformation to apply. Supported values:
    • "percentage": For data bounded between 0 and 100 (e.g., percentages, rates)
    • "positive": For strictly positive data (uses log transformation)
    • "boxcox": Applies Box-Cox transformation with automatically fitted λ parameter
  • values::Vector{F}: The input data values used to fit transformation parameters and determine offset

Returns

A tuple (forward_transform, inverse_transform) where:

  • forward_transform: Function that transforms data from original scale to transformed scale
  • inverse_transform: Function that transforms data from transformed scale back to original scale

Transformation Details

Percentage Transformation

  • Use case: Data bounded between 0 and 100 (percentages, rates)
  • Forward: y ↦ logit((y + offset) / 100)
  • Inverse: y ↦ max(logistic(y) * 100 - offset, 0)
  • Note: Uses logit/logistic to map [0,100] to (-∞,∞) and back

Positive Transformation

  • Use case: Strictly positive continuous data
  • Forward: y ↦ log(y + offset)
  • Inverse: y ↦ max(exp(y) - offset, 0)
  • Note: Log transformation for positive data with offset for numerical stability

Box-Cox Transformation

  • Use case: General purpose transformation for positive data
  • Forward: y ↦ BoxCox_λ(y + offset) where λ is automatically fitted
  • Inverse: Custom inverse function handling edge cases for numerical stability
  • Note: Automatically determines optimal λ parameter via maximum likelihood
  • Fallback: On near-constant data the Box-Cox MLE can return a pathological λ that collapses the transform to a constant (singular GP covariance, issue #51). When the transformed spread degenerates relative to a plain log, this falls back to the "positive" (log) transformation and emits a warning.

Offset Calculation

An offset is automatically calculated using _get_offet(values):

  • If minimum value is 0: offset = (minimum positive value) / 2
  • Otherwise: offset = 0
  • Purpose: Ensures numerical stability and handles boundary cases

Examples

# Percentage data (0-100 range)
values = [10.5, 25.3, 67.8, 89.2]
forward, inverse = get_transformations("percentage", values)
transformed = forward.(values)
recovered = inverse.(transformed)

# Strictly positive data
values = [1.2, 3.4, 8.9, 15.6]
forward, inverse = get_transformations("positive", values)

# General positive data with automatic Box-Cox fitting
values = [0.1, 0.5, 2.3, 5.7, 12.1]
forward, inverse = get_transformations("boxcox", values)

Throws

  • AssertionError: If transform_name is not one of the supported transformation types
  • AssertionError: Via _get_offet if values is empty or contains negative values

See Also

  • _get_offset: Calculates the offset value for numerical stability
  • _inv_boxcox: Handles inverse Box-Cox transformation with edge case handling
source
NowcastAutoGP.make_and_fit_modelMethod
function make_and_fit_model(
    data::TData; n_particles, smc_data_proportion,
    flat_threshold = 1.0e-3, config = AutoGP.GP.GPConfig(), kwargs...
)

Create and fit a Gaussian Process (GP) model using Sequential Monte Carlo (SMC) sampling.

Arguments

  • data: A data structure containing the dataset (data.ds) and the target values (data.y).
  • n_particles: The number of particles to use in the SMC sampling.
  • smc_data_proportion: The proportion of the data to use in each SMC step.
  • flat_threshold: If the (transformed) target's relative range is below this, small Gaussian jitter is added so the GP covariance stays positive-definite (default: 1.0e-3). Series with enough spread are left untouched, preserving existing behaviour.
  • config: An AutoGP.GP.GPConfig describing the GP prior — kernel-structure distribution (node_dist_leaf/node_dist_nocp/node_dist_cp) and hyperparameter priors (prior[:gamma]/prior[:period]/prior[:wildcard]). Defaults to GPConfig() (AutoGP defaults), so existing behaviour is unchanged. Build a customised one with the GPConfig keyword constructor — see Examples.
  • kwargs...: Additional keyword arguments forwarded to AutoGP.fit_smc!. Note that fit_smc! requires n_mcmc and n_hmc (number of MCMC structure steps and HMC parameter steps per SMC step), so these must be supplied here; other fit_smc! options (e.g. hmc_config, biased, shuffle, verbose) may also be passed.

Returns

  • model: The fitted GP model.

Examples

# Default prior:
model = make_and_fit_model(data; n_mcmc = 200, n_hmc = 50)

# Custom prior — pass a GPConfig. Kernel-structure fields go straight to the constructor,
# e.g. give SquaredExponential (code 3) leaf-kernel prior mass:
config = GPConfig(node_dist_leaf = [0.0, 0.25, 0.25, 0.25, 0.25])
model = make_and_fit_model(data; config = config, n_mcmc = 200, n_hmc = 50)

# To tweak a single nested prior entry (e.g. re-centre the periodic component near an annual
# cycle) copy the default prior and update the leaf:
prior = deepcopy(GPConfig().prior)
prior[:period][:mu] = log(1.0)
model = make_and_fit_model(data; config = GPConfig(prior = prior), n_mcmc = 200, n_hmc = 50)

@set from Accessors.jl is a handy alternative for these copy-and-update edits if you have it in your environment (config = @set GPConfig().prior[:period][:mu] = log(1.0)); it is not a dependency of this package.

source