derivkit.forecasting.laplace module#

Laplace approximation utilities.

Implements a Gaussian approximation to a posterior around a point theta_map. This is typically the maximum a posteriori (MAP). It uses a second-order Taylor expansion of the negative log-posterior:

neg_logpost(theta) = -logposterior(theta)
neg_logpost(theta) ≈ neg_logpost(theta_map) + 0.5 * d^T H d

where d = theta - theta_map and H is the Hessian of neg_logpost evaluated at theta_map. The approximate covariance is cov H^{-1}.

derivkit.forecasting.laplace.laplace_approximation(*, neg_logposterior: Callable[[ndarray[tuple[Any, ...], dtype[float64]]], float], theta_map: Sequence[float] | ndarray[tuple[Any, ...], dtype[float64]], method: str | None = None, n_workers: int = 1, dk_kwargs: Mapping[str, Any] | None = None, ensure_spd: bool = True, rcond: float = 1e-12) dict[str, Any]#

Computes a Laplace (Gaussian) approximation around theta_map.

The Laplace approximation replaces the posterior near its peak with a Gaussian. It does this by measuring the local curvature of the negative log-posterior using its Hessian at theta_map. The Hessian acts like a local precision matrix, and its inverse is the approximate covariance.

This is useful when a fast, local summary of the posterior is needed without running a full MCMC sampler. The output includes the expansion point, negative log-posterior value there, the Hessian, and the covariance matrix. In this approximation the Hessian acts as a local precision matrix, and its inverse is the approximate covariance.

Parameters:
  • neg_logposterior – Callable that accepts a 1D float64 parameter vector and returns a scalar negative log-posterior value.

  • theta_map – Expansion point for the approximation. This is often the maximum a posteriori estimate (MAP), which is the parameter vector that maximizes the posterior.

  • method – Derivative method name/alias forwarded to the Hessian builder.

  • n_workers – Outer parallelism forwarded to Hessian construction.

  • dk_kwargs – Extra keyword arguments forwarded to derivkit.DerivativeKit.differentiate().

  • ensure_spd – If True, attempt to regularize the Hessian to be symmetric positive definite (SPD) by adding diagonal jitter. This ensures the Gaussian approximation has a valid covariance matrix.

  • rcond – Cutoff for small singular values used by the pseudoinverse fallback when computing the covariance.

Returns:

  • “theta_map”: 1D array of the expansion point (dtype float64)).

  • ”neg_logposterior_at_map”: negative log-posterior at the expansion point (dtype float).

  • ”hessian”: array containing the Hessian of the negative log-posterior with shape (p, p) (local precision, dtype float64).

  • ”cov”: covariance matrix with shape (p, p) (approximate inverse Hessian).

  • ”jitter”: amount of diagonal jitter added (0.0 if None, dtype float).

Return type:

A dictionary with key-value pairs

Raises:
  • TypeError – If neg_logposterior is not scalar-valued.

  • ValueError – If inputs are invalid or non-finite values are encountered.

  • np.linalg.LinAlgError – If ensure_spd=True and the Hessian cannot be regularized to be SPD.

derivkit.forecasting.laplace.laplace_covariance(hessian: ndarray[tuple[Any, ...], dtype[float64]], *, rcond: float = 1e-12) ndarray[tuple[Any, ...], dtype[float64]]#

Computes the Laplace covariance matrix from a Hessian.

In the Laplace (Gaussian) approximation, the Hessian of the negative log-posterior at the expansion point acts like a local precision matrix. The approximate posterior covariance is the matrix inverse of that Hessian.

This helper inverts the Hessian using DerivKit’s robust linear-solve routine. It falls back to a pseudoinverse when the Hessian is ill-conditioned and returns a symmetrized covariance matrix suitable for downstream use, such as Gaussian approximations, uncertainty summaries, and proposal covariances.

Parameters:
  • hessian – 2D square Hessian matrix. This is usually the Hessian of the negative log-posterior evaluated at the expansion point.

  • rcond – Cutoff for small singular values used by the pseudoinverse fallback.

Returns:

A 2D symmetric covariance matrix with the same shape as hessian.

Raises:

ValueError – If hessian is not a finite square matrix.

derivkit.forecasting.laplace.laplace_hessian(*, neg_logposterior: Callable[[ndarray[tuple[Any, ...], dtype[float64]]], float], theta_map: Sequence[float] | ndarray[tuple[Any, ...], dtype[float64]], method: str | None = None, n_workers: int = 1, dk_kwargs: Mapping[str, Any] | None = None) ndarray[tuple[Any, ...], dtype[float64]]#

Computes the Hessian of the negative log-posterior at theta_map.

The Hessian at theta_map measures the local curvature of the posterior peak. In the Laplace approximation, this Hessian plays the role of a local precision matrix, and its inverse provides a fast Gaussian estimate of parameter uncertainties and correlations.

Internally, this uses derivkit.calculus_kit.CalculusKit (and therefore DerivKit’s Hessian construction machinery) to differentiate the scalar objective neg_logposterior(theta).

Parameters:
  • neg_logposterior – Callable returning the scalar negative log-posterior.

  • theta_map – Point where the curvature is evaluated (typically the MAP).

  • method – Derivative method name/alias forwarded to the calculus machinery.

  • n_workers – Outer parallelism forwarded to Hessian construction.

  • dk_kwargs – Extra keyword arguments forwarded to DerivativeKit.differentiate().

Returns:

A symmetric 2D array with shape (p, p) giving the Hessian of neg_logposterior evaluated at theta_map.

Raises:
  • TypeError – If neg_logposterior is not scalar-valued (Hessian is not 2D).

  • ValueError – If inputs are invalid or the Hessian is not a finite square matrix.

derivkit.forecasting.laplace.negative_logposterior(theta: Sequence[float] | ndarray[tuple[Any, ...], dtype[float64]], *, logposterior: Callable[[ndarray[tuple[Any, ...], dtype[float64]]], float]) float#

Computes the negative log-posterior at theta.

This converts a user-provided log-posterior (often returned as a log-likelihoods plus log-prior) into the objective most optimizers and curvature-based methods work with: the negative log-posterior.

In practice, this is the scalar function whose Hessian at the maximum a posteriori (MAP) point defines the Laplace (Gaussian) approximation to the posterior, and it is also the quantity minimized by MAP estimation routines.

Parameters:
  • theta – 1D array-like parameter vector.

  • logposterior – Callable that accepts a 1D float64 array and returns a scalar float.

Returns:

Negative log-posterior value as a float.

Raises:

ValueError – If theta is not a finite 1D vector or if the negative log-posterior evaluates to a non-finite value.