derivkit.forecast_kit module#
Provides the ForecastKit class.
A light wrapper around the core forecasting utilities
(derivkit.forecasting.fisher.build_fisher_matrix(),
derivkit.forecasting.dali.build_dali(),
derivkit.forecasting.fisher.build_delta_nu(),
and derivkit.forecasting.fisher.build_fisher_bias()) that exposes a simple
API for Fisher and DALI tensors.
Typical usage example:
>>> import numpy as np
>>> from derivkit.forecast_kit import ForecastKit
>>>
>>> # Toy linear model: 2 params -> 2 observables
>>> def model(theta: np.ndarray) -> np.ndarray:
... theta = np.asarray(theta, dtype=float)
... return np.array([theta[0] + 2.0 * theta[1], 3.0 * theta[0] - theta[1]], dtype=float)
>>>
>>> theta0 = np.array([0.1, -0.2])
>>> cov = np.eye(2)
>>>
>>> fk = ForecastKit(function=model, theta0=theta0, cov=cov)
>>> fisher_matrix = fk.fisher(method="finite", n_workers=1)
>>> fisher_matrix.shape
(2, 2)
>>>
>>> data_unbiased = model(theta0)
>>> data_biased = data_unbiased + np.array([1e-3, -2e-3])
>>> dn = fk.delta_nu(data_biased=data_biased, data_unbiased=data_unbiased)
>>> dn.shape
(2,)
>>>
>>> bias_vec, delta_theta = fk.fisher_bias(
... fisher_matrix=fisher_matrix,
... delta_nu=dn,
... method="finite",
... n_workers=1,
... )
>>> bias_vec.shape, delta_theta.shape
((2,), (2,))
- class derivkit.forecast_kit.ForecastKit(function: Callable[[Sequence[float] | ndarray], ndarray] | None, theta0: Sequence[float] | ndarray, cov: ndarray | Callable[[ndarray], ndarray])#
Bases:
objectProvides access to Fisher and DALI likelihoods-expansion tensors.
Initialises the ForecastKit with model, fiducials, and covariance.
- Parameters:
function – Callable returning the model mean vector \(\mu(\theta)\). May be
Noneif you only plan to use covariance-only workflows (e.g. generalized Fisher withterm="cov"). Required forfisher(),dali(), andfisher_bias().theta0 – Fiducial parameter values of shape
(p,)wherepis the number of parameters.cov –
Covariance specification. Supported forms are:
cov=C0: fixed covariance matrix \(C(\theta_0)\) with shape(n_obs, n_obs), wheren_obsis the number of observables.cov=cov_fn: callable withcov_fn(theta)returning the covariance matrix \(C(\theta)\) evaluated at the parameter vectortheta, with shape(n_obs, n_obs). The covariance attheta0is evaluated once and cached.
- dali(*, method: str | None = None, n_workers: int = 1, **dk_kwargs: Any) tuple[ndarray, ndarray]#
Builds the doublet-DALI tensors (G, H) for the given model.
- Parameters:
method – Method name or alias (e.g.,
"adaptive","finite"). IfNone, thederivkit.derivative_kit.DerivativeKitdefault is used.n_workers – Number of workers for per-parameter parallelization/threads. Default
1(serial). Inner batch evaluation is kept serial to avoid oversubscription.dk_kwargs – Additional keyword arguments passed to
derivkit.calculus_kit.CalculusKit.
- Returns:
A tuple
(G, H)whereGhas shape(p, p, p)andHhas shape(p, p, p, p), withpbeing the number of model parameters.
- delta_nu(data_unbiased: ndarray, data_biased: ndarray)#
Computes the difference between two data vectors.
This helper is used in Fisher-bias calculations and any other workflow where two data vectors are compared: it takes a pair of vectors (for example, a version with a systematic and one without) and returns their difference as a 1D array whose length matches the number of observables implied by
cov. It works with both 1D inputs and 2D arrays (for example, correlation-by-ell) and flattens 2D inputs using NumPy’s row-major (“C”) order, which is the standard convention throughout the DerivKit package.- Parameters:
data_unbiased – Reference data vector without the systematic. Can be 1D or 2D. If 1D, it must follow the NumPy’s row-major (“C”) flattening convention used throughout the package.
data_biased – Data vector that includes the systematic effect. Can be 1D or 2D. If 1D, it must follow the NumPy’s row-major (“C”) flattening convention used throughout the package.
- Returns:
A 1D NumPy array of length
n_observablesrepresenting the mismatch between the two input data vectors. This is simply the element-wise difference between the input with systematic and the input without systematic, flattened if necessary to match the expected observable ordering.- Raises:
ValueError – If input shapes differ, inputs are not 1D/2D, or the flattened length does not match
n_observables.FloatingPointError – If non-finite values are detected in the result.
- fisher(*, method: str | None = None, n_workers: int = 1, **dk_kwargs: Any) ndarray#
Computes the Fisher information matrix for a given model and covariance.
- Parameters:
method – Derivative method name or alias (e.g.,
"adaptive","finite"). IfNone, thederivkit.derivative_kit.DerivativeKitdefault is used.n_workers – Number of workers for per-parameter parallelisation. Default is
1(serial).**dk_kwargs – Additional keyword arguments forwarded to
derivkit.derivative_kit.DerivativeKit.differentiate().
- Returns:
Fisher matrix with shape
(n_parameters, n_parameters).
- fisher_bias(*, fisher_matrix: ndarray, delta_nu: ndarray, method: str | None = None, n_workers: int = 1, rcond: float = 1e-12, **dk_kwargs: Any) tuple[ndarray, ndarray]#
Estimates parameter bias using the stored model, expansion point, and covariance.
This function takes a model, an expansion point, a covariance matrix, a Fisher matrix, and a data-vector difference
delta_nuand maps that difference into parameter space. A common use case is the classic “Fisher bias” setup, where one asks how a systematic-induced change in the data would shift inferred parameters.Internally, the function evaluates the model response at the expansion point and uses the covariance and Fisher matrix to compute both the parameter-space bias vector and the corresponding shifts. See https://arxiv.org/abs/0710.5171 for details.
- Parameters:
fisher_matrix – Square matrix describing information about the parameters. Its shape must be
(p, p), wherepis the number of parameters.delta_nu – Difference between a biased and an unbiased data vector, for example \(\Delta\nu = \nu_{\mathrm{biased}} - \nu_{\mathrm{unbiased}}\). Accepts a 1D array of length n or a 2D array that will be flattened in row-major order (“C”) to length n, where n is the number of observables. If supplied as a 1D array, it must already follow the same row-major (“C”) flattening convention used throughout the package.
n_workers – Number of workers used by the internal derivative routine when forming the Jacobian.
method – Method name or alias (e.g.,
"adaptive","finite"). IfNone, thederivkit.derivative_kit.DerivativeKitdefault is used.rcond – Regularization cutoff for pseudoinverse. Default is
1e-12.**dk_kwargs – Additional keyword arguments passed to
derivkit.derivative_kit.DerivativeKit.differentiate().
- Returns:
A tuple
(bias_vec, delta_theta)of 1D arrays with lengthp, wherebias_vecis the parameter-space bias vector anddelta_thetaare the corresponding parameter shifts.- Raises:
ValueError – If input shapes are inconsistent with the stored model, covariance, or the Fisher matrix dimensions.
FloatingPointError – If the difference vector contains
NaN.
- gaussian_fisher(*, method: str | None = None, n_workers: int = 1, rcond: float = 1e-12, symmetrize_dcov: bool = True, **dk_kwargs: Any) ndarray#
Computes the generalized Fisher matrix for parameter-dependent mean and/or covariance.
This function computes the generalized Fisher matrix for a Gaussian likelihoods with parameter-dependent mean and/or covariance. Uses
derivkit.forecasting.fisher_general.build_generalized_fisher_matrix().- Parameters:
method – Derivative method name or alias (e.g.,
"adaptive","finite").n_workers – Number of workers for per-parameter parallelisation.
rcond – Regularization cutoff for pseudoinverse fallback in linear solves.
symmetrize_dcov – If
True, symmetrize each covariance derivative via \(\tfrac{1}{2}(C_{,i} + C_{,i}^{\mathsf{T}})\).**dk_kwargs – Forwarded to the internal derivative calls.
- Returns:
Fisher matrix with shape
(p, p).