psiphy.forecasting

psiphy.forecasting.fisher.tqdm(iterable, **kwargs)[source]
class psiphy.forecasting.fisher.FisherMatrix(simulator, theta_fid, param_names, X_mean=0.0, X_std=1.0, verbose=True, cache_file=None, **sim_kwargs)[source]

Numerical Fisher matrix forecast.

F_ij = (dmu/dtheta_i)^T C^{-1} (dmu/dtheta_j)

Parameters:
  • simulator (callable) – Function with signature simulator(theta, initial_seed, noise_seed, **sim_kwargs) -> np.ndarray. Must return a 1D array (the summary statistic).

  • theta_fid (array-like, shape (n_params,)) – Fiducial parameter values.

  • param_names (list of str) – Names for each parameter.

  • X_mean (float or np.ndarray) – Mean for normalisation (0.0 = no normalisation).

  • X_std (float or np.ndarray) – Std for normalisation (1.0 = no normalisation).

  • verbose (bool) – Print progress messages.

  • cache_file (str or None) – Path to pickle file for caching simulation outputs across runs.

  • **sim_kwargs – Additional keyword arguments forwarded to the simulator.

Examples

>>> fm = FisherMatrix(my_sim, theta_fid=[0.0, 1.0], param_names=['mu', 'sigma'])
>>> fm.estimate_covariance(seeds=range(1, 11), noise_seeds=range(10))
>>> fm.estimate_derivatives(delta_frac=0.01)
>>> fm.compute()
>>> fm.plot_ellipses()
run_simulations(theta=None, seeds=range(1, 11), noise_seeds=range(0, 10), save=True)[source]

Pre-run and cache simulations.

Parameters:
  • theta (array-like or None) – Parameter values. Defaults to theta_fid.

  • seeds (iterable) – Initial condition seeds.

  • noise_seeds (iterable) – Noise seeds.

  • save (bool) – Save cache to file after running.

Return type:

self

cache_info()[source]

Print a summary of cached simulations grouped by theta.

test_simulator(seed=1, noise_seed=0)[source]

Run the simulator once at fiducial to verify it works.

Returns:

Normalised output at fiducial parameters.

Return type:

np.ndarray

estimate_covariance(seeds=range(1, 11), noise_seeds=range(0, 10), n_jobs=1)[source]

Estimate the data covariance matrix at the fiducial point.

Realisations are appended to any existing ones, so you can call this multiple times with different seed ranges to increase the sample size.

Parameters:
  • seeds (iterable) – Initial condition / cosmic variance seeds.

  • noise_seeds (iterable) – Instrument noise seeds.

  • n_jobs (int) – Number of parallel workers (joblib). -1 uses all CPUs. Each worker handles all noise_seeds for one IC seed, so the simulator’s in-memory coeval cache is reused within each job.

Return type:

self

add_modelling_error(mod_err=0.1)[source]

Add a fractional modelling error to the covariance.

Appends C_mod = diag(mod_err * mu_fid)^2 to the existing covariance and recomputes C_inv with the Hartlap correction.

Parameters:

mod_err (float or np.ndarray) – Fractional error per data bin.

Return type:

self

estimate_derivatives(delta_frac=0.01, seeds=range(1, 11), noise_seed=0, n_stencil=2, use_crn=False, n_jobs=1)[source]

Estimate dmu/dtheta via central finite differences.

Parameters:
  • delta_frac (float) – Fractional step size h relative to the fiducial value.

  • seeds (iterable) – Seeds to average over for stable derivatives.

  • noise_seed (int or None) – Noise seed for derivatives (0 = noiseless; None = average over 10).

  • n_stencil (int) – Stencil order: 2, 4, or 6. Higher order reduces truncation error O(h^n).

  • use_crn (bool) – Use Common Random Numbers to cancel cosmic variance in the subtraction.

  • n_jobs (int) – Number of parallel workers (joblib). -1 uses all CPUs. Each worker handles all noise_seeds for one (theta_step, IC seed) pair.

Return type:

self

convergence_test_delta_frac(delta_fracs=None, seeds=range(1, 11), noise_seed=0, n_stencil=2, use_crn=False, plot=True, restore=True)[source]

Test sensitivity of dmu/dtheta to the finite-difference step size.

Sweeps delta_fracs and records the per-parameter norm of the derivative and the 1D marginal constraint. A good step size shows a flat plateau in the middle: upturn at small h (stochastic noise) and at large h (truncation error).

Parameters:
  • delta_fracs (array-like or None) – Step sizes to test. Default: 10 log-spaced values from 1e-3 to 0.3.

  • seeds (iterable) – Forwarded to estimate_derivatives.

  • noise_seed (int or None) – Forwarded to estimate_derivatives.

  • n_stencil (int) – Forwarded to estimate_derivatives.

  • use_crn (bool) – Forwarded to estimate_derivatives.

  • plot (bool) – If True, show diagnostic plots.

  • restore (bool) – Restore self.dmu_dtheta to its value before the test.

Returns:

Per-parameter dict with keys ‘delta_fracs’, ‘norms’, ‘sigma_1d’.

Return type:

dict

compute(method='standard', n_bootstrap=0, n_split=None, shrinkage_target='diagonal')[source]

Compute the Fisher matrix. Results stored in self.F and self.F_inv.

Parameters:
  • method (str) – ‘standard’ — F = J C⁻¹ Jᵀ (default). ‘lfim’ — Score-trick estimator via Synthetic Likelihood. ‘CoultonWandelt2023’ — Geometric combined estimator (Eq. 18). ‘Shrinkage’ — Ledoit-Wolf-style covariance regularisation.

  • n_bootstrap (int) – (lfim only) Bootstrap resamples for uncertainty; populates self.F_std.

  • n_split (int or None) – (CoultonWandelt2023 only) Size of set A. Defaults to 50/50 split.

  • shrinkage_target (str) – (Shrinkage only) ‘diagonal’ or ‘identity’.

Return type:

self

sigma_1d()[source]

1D marginal 1-sigma errors (marginalised over all other parameters).

Returns:

sigma_i = sqrt(F_inv[i,i])

Return type:

np.ndarray, shape (n_params,)

sigma_2d()[source]

2D marginal errors for all parameter pairs.

Return type:

dict keyed by (i, j) → (sigma_i, sigma_j)

sigma_conditional()[source]

Conditional errors (all other parameters perfectly known).

Returns:

sigma_i = 1 / sqrt(F[i,i])

Return type:

np.ndarray, shape (n_params,)

correlation_matrix()[source]

Parameter correlation matrix derived from F_inv.

Return type:

np.ndarray, shape (n_params, n_params)

sub_fisher(param_indices)[source]

Extract a Fisher sub-matrix for a subset of parameters.

Parameters:

param_indices (list of int)

Return type:

F_sub, F_sub_inv, sigma_sub

ellipse_params(i, j, n_sigma=1)[source]

Compute ellipse parameters for the 2D marginal of parameters (i, j).

Parameters:
  • i (int) – Parameter indices.

  • j (int) – Parameter indices.

  • n_sigma (int) – 1 or 2.

Returns:

dict with keys

Return type:

center, width, height, angle_deg

plot_ellipses(others=None, labels=None, label_names=None, prior_range=None, n_sigma_list=(1, 2), figsize=None, colors=None, filled=False)[source]

Triangle plot of Fisher confidence ellipses.

Parameters:
  • others (dict or None) – Additional FisherMatrix objects to overlay, keyed by label string.

  • labels (str or None) – Label for this Fisher matrix in the legend.

  • label_names (dict or None) – param_name -> LaTeX label string (without $).

  • prior_range (dict or None) – param_name -> (low, high) for axis limits.

  • n_sigma_list (list of int) – Contour levels (1 and/or 2).

  • figsize (tuple or None) – Figure size.

  • colors (dict or None) – label -> colour.

  • filled (bool) – Fill the innermost ellipse.

Return type:

matplotlib.figure.Figure

plot_derivatives(label_names=None, figsize=None)[source]

Plot dmu/dtheta for each parameter.

Return type:

matplotlib.figure.Figure

plot_covariance(figsize=(6, 5))[source]

Plot the covariance and correlation matrices.

Return type:

matplotlib.figure.Figure

sample(n_samples=100000)[source]

Draw Gaussian samples from the Fisher posterior.

Parameters:

n_samples (int)

Return type:

np.ndarray, shape (n_samples, n_params)

to_getdist(n_samples=100000, label=None, label_names=None)[source]

Generate a getdist MCSamples object (Gaussian approximation).

Parameters:
  • n_samples (int)

  • label (str or None)

  • label_names (dict or None) – param_name -> LaTeX label.

Return type:

getdist.MCSamples

save(filepath)[source]

Save Fisher matrix results to an npz file.

Parameters:

filepath (str)

classmethod load(filepath, simulator=None, **sim_kwargs)[source]

Load a FisherMatrix from an npz file.

Parameters:
  • filepath (str)

  • simulator (callable or None)

Return type:

FisherMatrix