API Reference
Public API
TemporalDisaggregations.disaggregate Function
disaggregate(method, aggregate_values, interval_start, interval_end;
loss_norm=L2DistLoss(), output_period=Month(1), output_start=nothing,
output_end=nothing, weights=nothing, irls_tol=1e-8, irls_max_iter=50)Reconstruct an instantaneous time series from interval-averaged observations.
Arguments
method::DisaggregationMethod: Algorithm configuration. One of:Spline(; smoothness=1e-1, penalty_order, tension)Sinusoid(; smoothness_interannual)GP(; kernel, obs_noise, n_quad)
aggregate_values: Vector of n observed averages over each interval.interval_start,interval_end: Interval boundaries asDate/DateTime.loss_norm::DistanceLoss = L2DistLoss(): Loss function from LossFunctions.jl. Common choices:L2DistLoss()(least squares),L1DistLoss()(robust to outliers),HuberLoss(δ)(hybrid - L2 for small residuals, L1 for large). For Huber loss, specify the threshold δ (default 1.345 achieves 95% efficiency at the normal distribution). Note: Robust losses (L1, Huber) use IRLS with a fixed penalty term. The samesmoothnessparameter may produce different effective smoothness for different loss functions. Tunesmoothnessseparately for each loss type if matching visual smoothness is important.output_period::Dates.Period = Month(1): Output grid spacing.output_start: Grid anchorDateorDateTime(defaultnothing).output_end: Last date of the output grid asDateorDateTime. Defaults to the end of the data domain.weights: Optional vector of n positive per-observation weights (e.g.1 ./ σ²).nothing(default) uses uniform weights. For robust losses (L1, Huber), these are multiplied element-wise with the IRLS weights at each iteration.irls_tol::Float64 = 1e-8: Convergence tolerance for IRLS iterations (only used for non-L2 losses). Smaller values (e.g., 1e-10) give more accurate solutions but take longer; larger values (e.g., 1e-6) converge faster but may be less accurate. The tolerance controls the relative parameter change: iterations stop whenmax(|x_new - x|) / (‖x‖ + 1e-10) < irls_tol.irls_max_iter::Int = 50: Maximum number of IRLS iterations (only used for non-L2 losses). Increase if IRLS does not converge within 50 iterations for difficult problems; decrease for faster (but potentially less accurate) solutions.
Returns
DimStack with :signal and :std layers indexed by Ti(dates). For all methods, :std is the spatially-varying sandwich standard deviation: std(t*) = σ̂ · sqrt(q(t*)), where σ̂ is the weighted residual RMS of predicted vs. observed interval averages and q(t*) is a dimensionless coverage factor that is smaller where observations are dense and larger where they are sparse.
Examples
using LossFunctions
result = disaggregate(Spline(smoothness=1e-3), y, t1, t2)
result = disaggregate(GP(obs_noise=4.0), y, t1, t2; loss_norm=L1DistLoss())
result = disaggregate(Sinusoid(), y, t1, t2; output_period=Day(1))
result = disaggregate(Spline(), y, t1, t2; output_end=Date(2020, 6, 1))
result = disaggregate(Spline(), y, t1, t2; weights = 1 ./ σ²_obs)
result = disaggregate(Spline(), y, t1, t2; loss_norm=HuberLoss(1.345))
result = disaggregate(Spline(), y, t1, t2; loss_norm=HuberLoss(2.0))
# Fast L1 solution with looser tolerance (2-3× faster)
result = disaggregate(Spline(), y, t1, t2; loss_norm=L1DistLoss(), irls_tol=1e-6)
# High-precision L1 solution with tighter tolerance
result = disaggregate(Spline(), y, t1, t2; loss_norm=L1DistLoss(), irls_tol=1e-10)
# Fast L1 with fewer IRLS iterations (for very large problems)
result = disaggregate(Spline(), y, t1, t2; loss_norm=L1DistLoss(), irls_max_iter=20)
# More iterations for difficult convergence cases
result = disaggregate(Spline(), y, t1, t2; loss_norm=HuberLoss(1.0), irls_max_iter=100)Missing docstring.
Missing docstring for yeardecimal. Check Documenter's build log for details.
TemporalDisaggregations.DisaggregationMethod Type
DisaggregationMethodAbstract supertype for temporal disaggregation algorithms. Subtypes hold algorithm-specific parameters with sensible defaults.
sourceTemporalDisaggregations.Spline Type
Spline(; smoothness=1e-3, penalty_order=3, tension=0.0)Quartic B-spline antiderivative fit with P-spline regularization.
Models the antiderivative F(t) as a B-spline so that F(t2ᵢ) − F(t1ᵢ) equals the observed area for each interval; the instantaneous signal is x(t) = F′(t).
Knots are placed at fixed monthly spacing (~12 per year) to maintain consistent smoothness behavior. The spline is evaluated at the user-requested output_period, which can be finer (e.g., daily) or coarser (e.g., quarterly) than the knot spacing.
The :std layer in the returned DimStack is a spatially-varying sandwich standard deviation: lower where observations are dense, higher where they are sparse.
Keywords
smoothness::Float64 = 1e-1: Regularization strength λ (larger = smoother). Note: When using robust losses likeL1DistLoss(), the samesmoothnessvalue may produce smoother results than withL2DistLoss(). This is expected behavior — tunesmoothnessseparately for each loss function.penalty_order::Int = 3: Order of the difference penalty.tension::Float64 = 0.0: Tension penalty strength (0 = standard P-spline).
TemporalDisaggregations.Sinusoid Type
Sinusoid(; smoothness_interannual=1e-2)Parametric model: mean + trend + per-year anomalies + annual sinusoid. All interval integrals are solved analytically (no quadrature) — the fastest method.
The :std layer in the returned DimStack is a spatially-varying sandwich standard deviation: lower where observations are dense, higher where they are sparse.
Keywords
smoothness_interannual::Float64 = 1e-2: Ridge penalty on inter-annual anomalies γ.
TemporalDisaggregations.GP Type
GP(; kernel=with_lengthscale(Matern52Kernel(), 1/6), obs_noise=1.0, n_quad=5)Sparse inducing-point Gaussian Process (DTC approximation) on a monthly grid. Supports arbitrary KernelFunctions.jl kernels (sums, products, periodic, etc.).
The :std layer in the returned DimStack is a spatially-varying sandwich standard deviation: lower where observations are dense, higher where they are sparse.
Keywords
kernel: AnyKernelFunctions.jlkernel. Default: Matérn-5/2 with 2-month lengthscale.obs_noise::Float64 = 1.0: Observation noise variance σ² in the same units as y². Controls the GP posterior smoothness: smaller values → tighter fit to observations.n_quad::Int = 5: Gauss-Legendre quadrature points per interval.
TemporalDisaggregations.interval_average Function
interval_average(result, t1, t2) → Vector{Float64}Compute the time-average of a disaggregated signal over each observation interval [t1[i], t2[i]] using trapezoidal integration on the high-resolution output grid.
Arguments
result: return value ofdisaggregate(has.signalfield with a:Tidimension)t1,t2: vectors of interval start/end times (DateTime or any type accepted byDateFormats.yeardecimal)
Returns
Vector{Float64} of length length(t1), each entry being the mean signal value over the corresponding interval. Intervals that fall entirely outside the output grid are filled with the nearest boundary value.
Author
Alex S. Gardner, JPL, Caltech.
source