calab¶
CaLab: calcium imaging analysis tools — deconvolution and data preparation.
Attributes¶
Classes¶
Result from bi-exponential kernel fitting. |
|
Full result from CaDecon (automated deconvolution via InDeCa algorithm). |
|
Ground truth for a single simulated cell. |
|
Configuration for automated CaDecon deconvolution via the bridge. |
|
Full result from FISTA deconvolution. |
|
Manages a Playwright browser for headless CaDecon / CaTune runs. |
|
Double-exponential kernel: h(t) = exp(-t/tau_decay) - exp(-t/tau_rise). |
|
Two-state HMM spike generator (silent/active) with bursty firing. |
|
Noise model: Gaussian + optional Poisson (shot) noise. |
|
Exponential photobleaching: |
|
Homogeneous Poisson spike generator. |
|
Mean-reverting Gaussian random walk baseline drift (default). |
|
Hill equation indicator saturation: F_sat = F^n / (F^n + Kd^n). |
|
Complete configuration for synthetic calcium trace generation. |
|
Complete simulation result with observed traces and per-cell ground truth. |
|
Deterministic sinusoidal baseline drift. |
|
Result from a single-trace InDeCa solve. |
|
Built-in indicator presets. Each method returns a SimulationConfig. |
Functions¶
|
Apply FFT bandpass filter derived from kernel time constants. Delegates to Rust. |
|
Build double-exponential calcium kernel. Delegates to Rust. |
|
Compute Lipschitz constant. Delegates to Rust. |
|
Compute the upsample factor for a given sampling rate and target. Delegates to Rust. |
|
Open CaDecon in the browser for automated deconvolution. |
|
Run deconvolution using parameters from a CaTune export JSON. |
|
Estimate a free-form kernel from traces and spike trains. Delegates to Rust. |
|
Fit a bi-exponential model to a free-form kernel. Delegates to Rust. |
|
Load traces from a CaImAn HDF5 results file. |
|
Load deconvolution parameters from a CaTune export JSON. |
|
Load traces from a Minian Zarr output directory. |
|
Load calcium traces and metadata saved by |
|
Run FISTA deconvolution on one or more calcium traces. |
|
Run FISTA deconvolution returning full results. |
|
Save calcium traces in CaTune-compatible format. |
|
Generate synthetic calcium imaging traces with full ground truth. |
|
Run the InDeCa pipeline on a single trace. Delegates to Rust. |
|
Derive AR(2) coefficients from tau parameters. |
|
Open CaTune in the browser for interactive parameter tuning. |
Package Contents¶
- class calab.BiexpFitResult¶
Bases:
NamedTupleResult from bi-exponential kernel fitting.
- class calab.CaDeconResult¶
Bases:
NamedTupleFull result from CaDecon (automated deconvolution via InDeCa algorithm).
- activity¶
Deconvolved activity matrix, shape
(n_cells, n_timepoints), float32.- Type:
np.ndarray
- alphas¶
Per-cell scaling factors, shape
(n_cells,), float64.- Type:
np.ndarray
- baselines¶
Per-cell baseline estimates, shape
(n_cells,), float64.- Type:
np.ndarray
- pves¶
Per-cell proportion of variance explained, shape
(n_cells,), float64.- Type:
np.ndarray
- kernel_slow¶
Slow biexponential kernel waveform, float32.
- Type:
np.ndarray
- kernel_fast¶
Fast biexponential kernel waveform, float32 (empty if single-component).
- Type:
np.ndarray
- activity: numpy.ndarray¶
- alphas: numpy.ndarray¶
- baselines: numpy.ndarray¶
- kernel_fast: numpy.ndarray¶
- kernel_slow: numpy.ndarray¶
- pves: numpy.ndarray¶
- class calab.CellGroundTruth¶
Bases:
pydantic.BaseModelGround truth for a single simulated cell.
- clean_calcium: numpy.ndarray¶
- model_config¶
- spikes: numpy.ndarray¶
- class calab.DeconConfig¶
Bases:
pydantic.BaseModelConfiguration for automated CaDecon deconvolution via the bridge.
Fields map 1:1 to the BridgeConfig TypeScript interface in @calab/io. Only
autorunis always serialized; all other fields useexclude_none=Trueso absent values fall through to browser defaults.
- class calab.DeconvolutionResult¶
Bases:
NamedTupleFull result from FISTA deconvolution.
- activity¶
Non-negative deconvolved activity estimates, same shape as input traces.
- Type:
np.ndarray
- reconvolution¶
K*activity + baseline, the model fit to the trace.
- Type:
np.ndarray
- converged¶
Whether convergence criterion was met (per-trace if multi-trace input).
- Type:
bool | np.ndarray
- activity: numpy.ndarray¶
- baseline: float | numpy.ndarray¶
- converged: bool | numpy.ndarray¶
- iterations: int | numpy.ndarray¶
- reconvolution: numpy.ndarray¶
- class calab.HeadlessBrowser(*, visible: bool = False)¶
Manages a Playwright browser for headless CaDecon / CaTune runs.
Can be used as a context manager for single runs, or kept alive across multiple
decon()/tune()calls for batch processing.Examples
Single run:
with HeadlessBrowser() as hb: result = calab.decon(traces, fs, headless=hb, autorun=True)
Batch (reuses browser across datasets):
with HeadlessBrowser() as hb: for traces in datasets: result = calab.decon(traces, fs, headless=hb, autorun=True)
- close() None¶
Shut down page, context, browser, and Playwright.
Each teardown runs in its own try/except so a crashed browser (where
context.close()raises) still letsbrowser.close()andpw.stop()run. Without this, a single teardown error would leak the remaining resources — most visibly an orphaned Chromium process that lingers until the Python process exits.
Navigate the managed page to url.
If the page is already on a different URL, it navigates there (effectively reusing the same tab).
- start() None¶
Launch the browser. Called automatically by
__enter__.If any initialization step raises after earlier steps succeeded,
close()runs before re-raising. This prevents partial-init failures (e.g. Chromium launch fails after the Playwright driver started) from leaking driver processes and Chromium instances across retries.
- property page: playwright.sync_api.Page¶
The managed Playwright
Page(raises if not started).
- class calab.KernelConfig¶
Bases:
pydantic.BaseModelDouble-exponential kernel: h(t) = exp(-t/tau_decay) - exp(-t/tau_rise).
Attribution: standard calcium response model (CaImAn, Suite2p, CaLab).
- class calab.MarkovConfig¶
Bases:
pydantic.BaseModelTwo-state HMM spike generator (silent/active) with bursty firing.
Attribution: CaLab web simulator Markov spike model.
- model_type: Literal['markov'] = 'markov'¶
- class calab.NoiseConfig¶
Bases:
pydantic.BaseModelNoise model: Gaussian + optional Poisson (shot) noise.
Attribution: Gaussian from CaLab web simulator. Shot noise from CASCADE (Rupprecht et al., 2021).
- class calab.PhotobleachingConfig¶
Bases:
pydantic.BaseModelExponential photobleaching:
F(t) *= 1 - amp * (1 - exp(-t/tau)).Attribution: NAOMi (Charles et al., 2019).
- class calab.PoissonConfig¶
Bases:
pydantic.BaseModelHomogeneous Poisson spike generator.
Attribution: standard model in OASIS (Friedrich et al., 2017) and CaImAn (Giovannucci et al., 2019).
- model_type: Literal['poisson'] = 'poisson'¶
- class calab.RandomWalkDrift¶
Bases:
pydantic.BaseModelMean-reverting Gaussian random walk baseline drift (default).
Models slow irregular baseline fluctuations from tissue movement, focus drift, and neuropil signal changes. From MLspike (Deneux et al., 2016, Nature Communications).
- model_type: Literal['random_walk'] = 'random_walk'¶
- class calab.SaturationConfig¶
Bases:
pydantic.BaseModelHill equation indicator saturation: F_sat = F^n / (F^n + Kd^n).
Attribution: MLspike (Deneux et al., 2016).
- class calab.SimulationConfig¶
Bases:
pydantic.BaseModelComplete configuration for synthetic calcium trace generation.
Per-cell variation (_cv fields) live on each config struct alongside the nominal value they modify. Alpha is here because it doesn’t belong to any pipeline step.
- drift: DriftModel¶
- kernel: KernelConfig¶
- noise: NoiseConfig¶
- photobleaching: PhotobleachingConfig¶
- saturation: SaturationConfig¶
- spike_model: SpikeModel¶
- class calab.SimulationResult¶
Bases:
pydantic.BaseModelComplete simulation result with observed traces and per-cell ground truth.
- config: SimulationConfig¶
- ground_truth: list[CellGroundTruth]¶
- model_config¶
- traces: numpy.ndarray¶
- class calab.SinusoidalDrift¶
Bases:
pydantic.BaseModelDeterministic sinusoidal baseline drift.
Useful as a simple test signal but not physically motivated.
- model_type: Literal['sinusoidal'] = 'sinusoidal'¶
- class calab.SolveTraceResult¶
Bases:
NamedTupleResult from a single-trace InDeCa solve.
- s_counts¶
Spike counts at the original sampling rate, shape
(n_timepoints,), float32.- Type:
np.ndarray
- s_counts: numpy.ndarray¶
- class calab.presets¶
Built-in indicator presets. Each method returns a SimulationConfig.
- static clean(**overrides: object) SimulationConfig¶
Near-ideal traces: minimal noise, no drift. For algorithm debugging.
- static gcamp6f(**overrides: object) SimulationConfig¶
GCaMP6f. Time constants from Chen et al., 2013, Nature.
- static gcamp6m(**overrides: object) SimulationConfig¶
GCaMP6m. Moderate kinetics. Chen et al., 2013.
- static gcamp6s(**overrides: object) SimulationConfig¶
GCaMP6s. Slow kinetics, high SNR. Chen et al., 2013.
- static jgcamp8f(**overrides: object) SimulationConfig¶
jGCaMP8f. Fast indicator, noisier. Zhang et al., 2023.
- static ogb1(**overrides: object) SimulationConfig¶
OGB-1 synthetic dye. Stosiek et al., 2003.
- calab.bandpass_filter(trace: numpy.ndarray, tau_rise: float, tau_decay: float, fs: float) numpy.ndarray¶
Apply FFT bandpass filter derived from kernel time constants. Delegates to Rust.
- calab.build_kernel(tau_rise: float, tau_decay: float, fs: float) numpy.ndarray¶
Build double-exponential calcium kernel. Delegates to Rust.
- calab.compute_lipschitz(kernel: numpy.ndarray) float¶
Compute Lipschitz constant. Delegates to Rust.
- calab.compute_upsample_factor(fs: float, target_fs: float) int¶
Compute the upsample factor for a given sampling rate and target. Delegates to Rust.
- calab.decon(traces: numpy.ndarray, fs: float = 30.0, timeout: float | None = None, port: int | None = None, app_url: str | None = None, open_browser: bool = True, headless: calab._bridge._headless.HeadlessBrowser | bool | None = None, *, autorun: bool = False, upsample_target: int | None = None, hp_filter_enabled: bool | None = None, lp_filter_enabled: bool | None = None, max_iterations: int | None = None, convergence_tol: float | None = None, num_subsets: int | None = None, target_coverage: float | None = None, aspect_ratio: float | None = None, seed: int | None = None) calab._compute.CaDeconResult | None¶
Open CaDecon in the browser for automated deconvolution.
Starts a localhost HTTP server serving the provided traces, opens CaDecon with a
?bridge=parameter pointing to the server, and waits for the browser to export deconvolution results back.- Parameters:
traces (np.ndarray) – Calcium traces, shape
(n_cells, n_timepoints)or(n_timepoints,).fs (float) – Sampling rate in Hz. Default: 30.0.
timeout (float, optional) – Seconds to wait for results. None = wait forever (until Ctrl-C).
port (int, optional) – Port to bind to. None = auto-assign.
app_url (str, optional) – Override CaDecon URL (for local dev). Default: GitHub Pages.
open_browser (bool) – Whether to auto-open the browser. Default: True.
headless (HeadlessBrowser or bool or None) –
None/False: default (usewebbrowser.open).True: create a temporary headless browser for this call.HeadlessBrowser: reuse an existing browser instance (for batch).autorun (bool) – If True, the solver starts automatically after loading. Default: False.
upsample_target (int, optional) – Target sampling rate for upsampling. Must be > 0.
hp_filter_enabled (bool, optional) – Enable high-pass filter.
lp_filter_enabled (bool, optional) – Enable low-pass filter.
max_iterations (int, optional) – Maximum solver iterations (1–200).
convergence_tol (float, optional) – Convergence tolerance (0–1 exclusive).
num_subsets (int, optional) – Number of random subsets. Must be > 0.
target_coverage (float, optional) – Target coverage fraction (0–1].
aspect_ratio (float, optional) – Subset aspect ratio. Must be > 0.
seed (int, optional) – Random seed for subset placement.
- Returns:
Deconvolution results if received, None if timeout/cancelled.
- Return type:
CaDeconResult or None
- calab.deconvolve_from_export(traces: numpy.ndarray, params_path: str | pathlib.Path, return_full: bool = False)¶
Run deconvolution using parameters from a CaTune export JSON.
Loads parameters via
load_export_params(), optionally applies the bandpass filter, and runs FISTA deconvolution.- Parameters:
traces (np.ndarray) – Calcium traces, shape
(n_timepoints,)or(n_cells, n_timepoints).params_path (str or Path) – Path to the CaTune export JSON file.
return_full (bool, optional) – If True, return a
DeconvolutionResultnamedtuple. Default False returns just the activity array.
- Returns:
Deconvolved activity (or full result if
return_full=True).- Return type:
np.ndarray or DeconvolutionResult
- calab.estimate_kernel(traces_flat: numpy.ndarray, spikes_flat: numpy.ndarray, trace_lengths: numpy.ndarray, alphas: numpy.ndarray, baselines: numpy.ndarray, kernel_length: int, *, max_iters: int = 200, tol: float = 0.0001, warm_kernel: numpy.ndarray | None = None, smooth_lambda: float = 0.0) numpy.ndarray¶
Estimate a free-form kernel from traces and spike trains. Delegates to Rust.
- Parameters:
traces_flat (np.ndarray) – Concatenated 1-D traces (all cells flattened).
spikes_flat (np.ndarray) – Concatenated 1-D spike trains (matching traces_flat).
trace_lengths (np.ndarray) – Length of each individual trace in the flat arrays.
alphas (np.ndarray) – Per-trace amplitude scaling factors.
baselines (np.ndarray) – Per-trace baseline estimates.
kernel_length (int) – Desired output kernel length in samples.
max_iters (int) – Maximum FISTA iterations for kernel estimation.
tol (float) – Convergence tolerance.
warm_kernel (np.ndarray, optional) – Kernel from a previous iteration for warm-start.
smooth_lambda (float) – Total-variation smoothness penalty weight.
- Returns:
Estimated free-form kernel, shape
(kernel_length,), float32.- Return type:
np.ndarray
- calab.fit_biexponential(h_free: numpy.ndarray, fs: float, *, refine: bool = True, skip: int = 0, warm: BiexpFitResult | None = None) BiexpFitResult¶
Fit a bi-exponential model to a free-form kernel. Delegates to Rust.
- Parameters:
h_free (np.ndarray) – Free-form kernel (1-D).
fs (float) – Sampling rate in Hz.
refine (bool) – Whether to refine with a fast (second) component.
skip (int) – Number of leading samples to skip in the fit.
warm (BiexpFitResult, optional) – Previous fit result for warm-start.
- Return type:
- calab.load_caiman(path: str, trace_key: str = 'estimates/C', fs_key: str = 'params/data/fr', fs: float | None = None) tuple[numpy.ndarray, dict]¶
Load traces from a CaImAn HDF5 results file.
- Parameters:
- Returns:
traces (np.ndarray) – Traces array, shape
(n_cells, n_timepoints), dtype float64.metadata (dict) – Metadata dict with keys:
source,sampling_rate_hz,num_cells,num_timepoints.
- Raises:
ImportError – If h5py is not installed.
FileNotFoundError – If the HDF5 file does not exist.
KeyError – If
trace_keyis not found in the file.
- calab.load_export_params(path: str | pathlib.Path) dict¶
Load deconvolution parameters from a CaTune export JSON.
Parses the JSON export schema 1.1.0 (defined in src/lib/export.ts) and returns the parameters needed for deconvolution.
- Parameters:
path (str or Path) – Path to the CaTune export JSON file.
- Returns:
Dictionary with keys:
tau_rise,tau_decay,lambda_,fs,filter_enabled.- Return type:
- Raises:
FileNotFoundError – If the JSON file does not exist.
KeyError – If required parameter fields are missing.
- calab.load_minian(path: str, trace_key: str = 'C', fs: float | None = None) tuple[numpy.ndarray, dict]¶
Load traces from a Minian Zarr output directory.
- Parameters:
- Returns:
traces (np.ndarray) – Traces array, shape
(n_cells, n_timepoints), dtype float64.metadata (dict) – Metadata dict with keys:
source,sampling_rate_hz,num_cells,num_timepoints.
- Raises:
ImportError – If zarr is not installed.
FileNotFoundError – If the Zarr directory does not exist.
KeyError – If
trace_keyis not found in the store.
- calab.load_tuning_data(path: str | pathlib.Path) tuple[numpy.ndarray, dict]¶
Load calcium traces and metadata saved by
save_for_tuning().- Parameters:
path (str or Path) – Path stem (without extension), matching the
pathargument used insave_for_tuning().- Returns:
traces (np.ndarray) – Loaded traces array, dtype float64.
metadata (dict) – Metadata from the JSON sidecar.
- Raises:
FileNotFoundError – If either
{path}.npyor{path}_metadata.jsonis missing.
- calab.run_deconvolution(traces: numpy.ndarray, fs: float, tau_r: float, tau_d: float, lam: float, max_iters: int = 2000, conv_mode: str = 'fft', constraint: str = 'nonneg') numpy.ndarray¶
Run FISTA deconvolution on one or more calcium traces.
Delegates to the Rust solver via calab._solver.
- Parameters:
traces (np.ndarray) – Input traces, shape
(n_timepoints,)for a single trace or(n_cells, n_timepoints)for multiple traces.fs (float) – Sampling rate in Hz.
tau_r (float) – Rise time constant in seconds.
tau_d (float) – Decay time constant in seconds.
lam (float) – L1 penalty (sparsity regularization strength).
max_iters (int, optional) – Maximum number of FISTA iterations, by default 2000.
conv_mode (str, optional) – Convolution mode:
'fft'(default) or'banded'(O(T) AR2).constraint (str, optional) – Constraint type:
'nonneg'(default, L1 + non-negative) or'box01'(box constraint [0, 1], no L1 penalty).
- Returns:
Non-negative activity estimates, same shape as input
traces.- Return type:
np.ndarray
- calab.run_deconvolution_full(traces: numpy.ndarray, fs: float, tau_r: float, tau_d: float, lam: float, max_iters: int = 2000, conv_mode: str = 'fft', constraint: str = 'nonneg') DeconvolutionResult¶
Run FISTA deconvolution returning full results.
- Parameters:
traces (np.ndarray) – Input traces, shape
(n_timepoints,)for a single trace or(n_cells, n_timepoints)for multiple traces.fs (float) – Sampling rate in Hz.
tau_r (float) – Rise time constant in seconds.
tau_d (float) – Decay time constant in seconds.
lam (float) – L1 penalty (sparsity regularization strength).
max_iters (int, optional) – Maximum number of FISTA iterations, by default 2000.
conv_mode (str, optional) – Convolution mode:
'fft'(default) or'banded'(O(T) AR2).constraint (str, optional) – Constraint type:
'nonneg'(default, L1 + non-negative) or'box01'(box constraint [0, 1], no L1 penalty).
- Returns:
Namedtuple with fields:
activity,baseline,reconvolution,iterations,converged.- Return type:
- calab.save_for_tuning(traces: numpy.ndarray, fs: float, path: str | pathlib.Path, metadata: dict | None = None) None¶
Save calcium traces in CaTune-compatible format.
- Creates two files:
{path}.npy– Float64 array, shape(n_cells, n_timepoints), C-contiguous, little-endian (<f8).{path}_metadata.json– JSON sidecar with sampling rate, schema version, dimensions, and optional user metadata.
The
.npyfile is loadable by CaTune’s browser.npyparser, which expectsdtype='<f8'andfortran_order=False.Converting from CaImAn (users have h5py):
import h5py import calab with h5py.File("caiman_results.hdf5", "r") as f: traces = f["estimates/C"][:] # shape: (n_cells, n_timepoints) fs = float(f["params/data/fr"][()]) calab.save_for_tuning(traces, fs, "my_recording") # -> my_recording.npy + my_recording_metadata.json
Converting from Minian (users have zarr/xarray):
import zarr import calab store = zarr.open("minian_output", mode="r") traces = store["C"][:] # shape: (n_cells, n_frames) fs = 30.0 # user must know their frame rate calab.save_for_tuning(traces, fs, "my_recording")
- Parameters:
traces (np.ndarray) – Calcium traces, shape
(n_timepoints,)for a single trace or(n_cells, n_timepoints)for multiple traces.fs (float) – Sampling rate in Hz.
path (str or Path) – Output path stem (without extension). Will create
{path}.npyand{path}_metadata.json.metadata (dict, optional) – Additional metadata to include in the JSON sidecar. Keys are merged into the output; built-in keys take precedence.
- Raises:
ValueError – If
traceshas more than 2 dimensions.
- calab.simulate(config: SimulationConfig | None = None, **kwargs: object) SimulationResult¶
Generate synthetic calcium imaging traces with full ground truth.
- Parameters:
config (SimulationConfig, optional) – Full configuration. If None, a default config is created.
**kwargs – Override fields on the default/provided config (e.g., num_cells=50, seed=123).
- Returns:
Contains traces array and per-cell ground truth.
- Return type:
- calab.solve_trace(trace: numpy.ndarray, tau_rise: float, tau_decay: float, fs: float, *, upsample_factor: int = 1, max_iters: int = 500, tol: float = 0.0001, hp_enabled: bool = False, lp_enabled: bool = False, warm_counts: numpy.ndarray | None = None, lambda_: float = 0.0) SolveTraceResult¶
Run the InDeCa pipeline on a single trace. Delegates to Rust.
- Parameters:
trace (np.ndarray) – 1-D calcium trace.
tau_rise (float) – Time constants in seconds.
tau_decay (float) – Time constants in seconds.
fs (float) – Sampling rate in Hz.
upsample_factor (int) – Upsampling multiplier (1 = no upsampling).
max_iters (int) – Maximum FISTA iterations.
tol (float) – Convergence tolerance.
hp_enabled (bool) – Enable high-pass / low-pass filtering.
lp_enabled (bool) – Enable high-pass / low-pass filtering.
warm_counts (np.ndarray, optional) – Spike counts from a previous iteration (at original rate) for warm-start.
lambda (float) – L1 sparsity penalty.
- Return type:
- calab.tau_to_ar2(tau_rise: float, tau_decay: float, fs: float) tuple[float, float, float, float]¶
Derive AR(2) coefficients from tau parameters.
Pure Python (trivial math, no solver needed).
- calab.tune(traces: numpy.ndarray, fs: float = 30.0, timeout: float | None = None, port: int | None = None, app_url: str | None = None, open_browser: bool = True) dict | None¶
Open CaTune in the browser for interactive parameter tuning.
Starts a localhost HTTP server serving the provided traces, opens CaTune with a
?bridge=parameter pointing to the server, and waits for the user to export parameters from the web app.- Parameters:
traces (np.ndarray) – Calcium traces, shape
(n_cells, n_timepoints)or(n_timepoints,).fs (float) – Sampling rate in Hz. Default: 30.0.
timeout (float, optional) – Seconds to wait for params. None = wait forever (until Ctrl-C).
port (int, optional) – Port to bind to. None = auto-assign.
app_url (str, optional) – Override CaTune URL (for local dev). Default: GitHub Pages.
open_browser (bool) – Whether to auto-open the browser. Default: True.
- Returns:
Exported parameters dict if received, None if timeout/cancelled. Keys:
tau_rise,tau_decay,lambda_,fs,filter_enabled.- Return type:
dict or None
- calab.DriftModel¶
- calab.SpikeModel¶