Skip to content

utils

Utility modules: quick-start presets, batched parameter sweeps, and HDF5 I/O.


Presets

Quick-start factory functions for common bubble configurations.

from jbubble.utils.presets import free_bubble, lipid_bubble, thick_shell_bubble

jbubble.utils.presets.BubblePreset

Bases: NamedTuple

Assembled (EoM, pulse) pair returned by preset factory functions.

Can be unpacked directly::

eom, pulse = free_bubble()

or accessed by name::

preset = free_bubble()
result = run_simulation(preset.eom, preset.pulse)

jbubble.utils.presets.free_bubble(R0=2e-06, freq=1000000.0, pressure=100000.0, cycle_num=5, *, gamma=1.4, sigma=0.072, mu=0.001, P_amb=101325.0, rho_L=998.0, c_L=1500.0)

Uncoated gas bubble in a Newtonian liquid (Keller-Miksis).

The simplest physically meaningful preset: a polytropic gas bubble with Laplace surface tension only — no shell coating. Suitable for free cavitation bubbles and as a baseline before adding a shell.

Physics: Keller-Miksis EoM + PolytropicGas + NoShell + NewtonianMedium.

Parameters:

Name Type Description Default
R0 float

Equilibrium radius [m]. Default: 2 µm.

2e-06
freq float

Driving frequency [Hz]. Default: 1 MHz.

1000000.0
pressure float

Peak acoustic pressure amplitude [Pa]. Default: 100 kPa.

100000.0
cycle_num int

Number of tone-burst cycles. Default: 5.

5
gamma float

Polytropic exponent. Default: 1.4 (diatomic gas / air).

1.4
sigma float

Surface tension [N/m]. Default: 0.072 (air–water interface).

0.072
mu float

Liquid dynamic viscosity [Pa·s]. Default: 1e-3 (water, 20 °C).

0.001
P_amb float

Ambient pressure [Pa]. Default: 101 325 (1 atm).

101325.0
rho_L float

Liquid density [kg/m³]. Default: 998 (water, 20 °C).

998.0
c_L float

Speed of sound in the liquid [m/s]. Default: 1500 (water).

1500.0

Returns:

Type Description
BubblePreset

(eom, pulse) ready for :func:~jbubble.run_simulation.

Source code in jbubble/utils/presets.py
def free_bubble(
    R0: float = 2e-6,
    freq: float = 1e6,
    pressure: float = 100e3,
    cycle_num: int = 5,
    *,
    gamma: float = 1.4,
    sigma: float = 0.072,
    mu: float = 1e-3,
    P_amb: float = 101325.0,
    rho_L: float = 998.0,
    c_L: float = 1500.0,
) -> BubblePreset:
    """Uncoated gas bubble in a Newtonian liquid (Keller-Miksis).

    The simplest physically meaningful preset: a polytropic gas bubble
    with Laplace surface tension only — no shell coating.  Suitable for
    free cavitation bubbles and as a baseline before adding a shell.

    Physics: Keller-Miksis EoM + PolytropicGas + NoShell + NewtonianMedium.

    Parameters
    ----------
    R0 : float
        Equilibrium radius [m].  Default: 2 µm.
    freq : float
        Driving frequency [Hz].  Default: 1 MHz.
    pressure : float
        Peak acoustic pressure amplitude [Pa].  Default: 100 kPa.
    cycle_num : int
        Number of tone-burst cycles.  Default: 5.
    gamma : float
        Polytropic exponent.  Default: 1.4 (diatomic gas / air).
    sigma : float
        Surface tension [N/m].  Default: 0.072 (air–water interface).
    mu : float
        Liquid dynamic viscosity [Pa·s].  Default: 1e-3 (water, 20 °C).
    P_amb : float
        Ambient pressure [Pa].  Default: 101 325 (1 atm).
    rho_L : float
        Liquid density [kg/m³].  Default: 998 (water, 20 °C).
    c_L : float
        Speed of sound in the liquid [m/s].  Default: 1500 (water).

    Returns
    -------
    BubblePreset
        ``(eom, pulse)`` ready for :func:`~jbubble.run_simulation`.
    """
    eom = KellerMiksis(
        gas=PolytropicGas(gamma=gamma),
        shell=NoShell(sigma=sigma),
        medium=NewtonianMedium(mu=mu),
        R0=R0,
        P_amb=P_amb,
        rho_L=rho_L,
        c_L=c_L,
    )
    pulse = ToneBurst(freq=freq, pressure=pressure, shape=Sine(), cycle_num=cycle_num)
    return BubblePreset(eom=eom, pulse=pulse)

jbubble.utils.presets.lipid_bubble(R0=2e-06, freq=1000000.0, pressure=100000.0, cycle_num=5, *, kappa_s=2.4e-09, chi=0.55, sigma_rupture=0.072, R_buckle_ratio=0.98, gamma=1.4, mu=0.001, P_amb=101325.0, rho_L=998.0, c_L=1500.0)

Lipid-shelled ultrasound contrast agent — Marmottant (2005) model.

Models a thin lipid monolayer with surface-dilatational viscosity and a smooth Gompertz surface tension law (preferred over the piecewise Marmottant law for gradient-based parameter fitting). Representative of clinical agents such as SonoVue and Definity.

Physics: Keller-Miksis EoM + PolytropicGas + LipidShell(GompertzSurfaceTension) + NewtonianMedium.

Parameters:

Name Type Description Default
R0 float

Equilibrium radius [m]. Default: 2 µm.

2e-06
freq float

Driving frequency [Hz]. Default: 1 MHz.

1000000.0
pressure float

Peak acoustic pressure amplitude [Pa]. Default: 100 kPa.

100000.0
cycle_num int

Number of tone-burst cycles. Default: 5.

5
kappa_s float

Shell surface-dilatational viscosity [N·s/m]. Default: 2.4e-9 (Marmottant 2005, BR14 / SonoVue-type).

2.4e-09
chi float

Shell elasticity [N/m]. Default: 0.55.

0.55
sigma_rupture float

Asymptotic (ruptured) surface tension [N/m]. Default: 0.072 (water).

0.072
R_buckle_ratio float

Buckling radius as fraction of R0. Default: 0.98 — gives an initial surface tension sigma(R0) ≈ 0.023 N/m, well below sigma_rupture, which is required for the Gompertz model to be well-posed.

0.98
gamma float

Polytropic exponent. Default: 1.4.

1.4
mu float

Liquid dynamic viscosity [Pa·s]. Default: 1e-3 (water).

0.001
P_amb float

Ambient pressure [Pa]. Default: 101 325 (1 atm).

101325.0
rho_L float

Liquid density [kg/m³]. Default: 998 (water).

998.0
c_L float

Speed of sound in the liquid [m/s]. Default: 1500 (water).

1500.0

Returns:

Type Description
BubblePreset

(eom, pulse) ready for :func:~jbubble.run_simulation.

References

Marmottant et al., J. Acoust. Soc. Am. 118 (2005) 3499–3505.

Source code in jbubble/utils/presets.py
def lipid_bubble(
    R0: float = 2e-6,
    freq: float = 1e6,
    pressure: float = 100e3,
    cycle_num: int = 5,
    *,
    kappa_s: float = 2.4e-9,
    chi: float = 0.55,
    sigma_rupture: float = 0.072,
    R_buckle_ratio: float = 0.98,
    gamma: float = 1.4,
    mu: float = 1e-3,
    P_amb: float = 101325.0,
    rho_L: float = 998.0,
    c_L: float = 1500.0,
) -> BubblePreset:
    """Lipid-shelled ultrasound contrast agent — Marmottant (2005) model.

    Models a thin lipid monolayer with surface-dilatational viscosity and a
    smooth Gompertz surface tension law (preferred over the piecewise
    Marmottant law for gradient-based parameter fitting).  Representative
    of clinical agents such as SonoVue and Definity.

    Physics: Keller-Miksis EoM + PolytropicGas
             + LipidShell(GompertzSurfaceTension) + NewtonianMedium.

    Parameters
    ----------
    R0 : float
        Equilibrium radius [m].  Default: 2 µm.
    freq : float
        Driving frequency [Hz].  Default: 1 MHz.
    pressure : float
        Peak acoustic pressure amplitude [Pa].  Default: 100 kPa.
    cycle_num : int
        Number of tone-burst cycles.  Default: 5.
    kappa_s : float
        Shell surface-dilatational viscosity [N·s/m].  Default: 2.4e-9
        (Marmottant 2005, BR14 / SonoVue-type).
    chi : float
        Shell elasticity [N/m].  Default: 0.55.
    sigma_rupture : float
        Asymptotic (ruptured) surface tension [N/m].  Default: 0.072 (water).
    R_buckle_ratio : float
        Buckling radius as fraction of R0.  Default: 0.98 — gives an
        initial surface tension ``sigma(R0) ≈ 0.023 N/m``, well below
        ``sigma_rupture``, which is required for the Gompertz model to
        be well-posed.
    gamma : float
        Polytropic exponent.  Default: 1.4.
    mu : float
        Liquid dynamic viscosity [Pa·s].  Default: 1e-3 (water).
    P_amb : float
        Ambient pressure [Pa].  Default: 101 325 (1 atm).
    rho_L : float
        Liquid density [kg/m³].  Default: 998 (water).
    c_L : float
        Speed of sound in the liquid [m/s].  Default: 1500 (water).

    Returns
    -------
    BubblePreset
        ``(eom, pulse)`` ready for :func:`~jbubble.run_simulation`.

    References
    ----------
    Marmottant et al., J. Acoust. Soc. Am. 118 (2005) 3499–3505.
    """
    shell = LipidShell(
        sigma=GompertzSurfaceTension(
            R_buckle_ratio=R_buckle_ratio,
            chi=chi,
            sigma_rupture=sigma_rupture,
        ),
        kappa_s=kappa_s,
    )
    eom = KellerMiksis(
        gas=PolytropicGas(gamma=gamma),
        shell=shell,
        medium=NewtonianMedium(mu=mu),
        R0=R0,
        P_amb=P_amb,
        rho_L=rho_L,
        c_L=c_L,
    )
    pulse = ToneBurst(freq=freq, pressure=pressure, shape=Sine(), cycle_num=cycle_num)
    return BubblePreset(eom=eom, pulse=pulse)

jbubble.utils.presets.thick_shell_bubble(R0=2e-06, freq=1000000.0, pressure=100000.0, cycle_num=5, *, d_s=1.5e-08, G_s=10000000.0, mu_s=0.5, sigma=0.04, gamma=1.4, mu=0.001, P_amb=101325.0, rho_L=998.0, c_L=1500.0)

Polymer thick-shell bubble — Church (1995) model.

Models a thick viscoelastic shell with elastic and viscous contributions. Representative of polymer-shelled agents such as Optison (albumin), or experimental PLGA microbubbles.

Physics: Keller-Miksis EoM + PolytropicGas + ThickShell + NewtonianMedium.

Parameters:

Name Type Description Default
R0 float

Equilibrium radius [m]. Default: 2 µm.

2e-06
freq float

Driving frequency [Hz]. Default: 1 MHz.

1000000.0
pressure float

Peak acoustic pressure amplitude [Pa]. Default: 100 kPa.

100000.0
cycle_num int

Number of tone-burst cycles. Default: 5.

5
d_s float

Shell thickness [m]. Default: 15 nm.

1.5e-08
G_s float

Shell shear modulus [Pa]. Default: 10 MPa (stiff polymer shell).

10000000.0
mu_s float

Shell viscosity [Pa·s]. Default: 0.5.

0.5
sigma float

Surface tension [N/m]. Default: 0.04 (reduced by polymer coating).

0.04
gamma float

Polytropic exponent. Default: 1.4.

1.4
mu float

Liquid dynamic viscosity [Pa·s]. Default: 1e-3 (water).

0.001
P_amb float

Ambient pressure [Pa]. Default: 101 325 (1 atm).

101325.0
rho_L float

Liquid density [kg/m³]. Default: 998 (water).

998.0
c_L float

Speed of sound in the liquid [m/s]. Default: 1500 (water).

1500.0

Returns:

Type Description
BubblePreset

(eom, pulse) ready for :func:~jbubble.run_simulation.

References

Church, J. Acoust. Soc. Am. 97 (1995) 1510–1521.

Source code in jbubble/utils/presets.py
def thick_shell_bubble(
    R0: float = 2e-6,
    freq: float = 1e6,
    pressure: float = 100e3,
    cycle_num: int = 5,
    *,
    d_s: float = 15e-9,
    G_s: float = 10e6,
    mu_s: float = 0.5,
    sigma: float = 0.04,
    gamma: float = 1.4,
    mu: float = 1e-3,
    P_amb: float = 101325.0,
    rho_L: float = 998.0,
    c_L: float = 1500.0,
) -> BubblePreset:
    """Polymer thick-shell bubble — Church (1995) model.

    Models a thick viscoelastic shell with elastic and viscous contributions.
    Representative of polymer-shelled agents such as Optison (albumin),
    or experimental PLGA microbubbles.

    Physics: Keller-Miksis EoM + PolytropicGas + ThickShell + NewtonianMedium.

    Parameters
    ----------
    R0 : float
        Equilibrium radius [m].  Default: 2 µm.
    freq : float
        Driving frequency [Hz].  Default: 1 MHz.
    pressure : float
        Peak acoustic pressure amplitude [Pa].  Default: 100 kPa.
    cycle_num : int
        Number of tone-burst cycles.  Default: 5.
    d_s : float
        Shell thickness [m].  Default: 15 nm.
    G_s : float
        Shell shear modulus [Pa].  Default: 10 MPa (stiff polymer shell).
    mu_s : float
        Shell viscosity [Pa·s].  Default: 0.5.
    sigma : float
        Surface tension [N/m].  Default: 0.04 (reduced by polymer coating).
    gamma : float
        Polytropic exponent.  Default: 1.4.
    mu : float
        Liquid dynamic viscosity [Pa·s].  Default: 1e-3 (water).
    P_amb : float
        Ambient pressure [Pa].  Default: 101 325 (1 atm).
    rho_L : float
        Liquid density [kg/m³].  Default: 998 (water).
    c_L : float
        Speed of sound in the liquid [m/s].  Default: 1500 (water).

    Returns
    -------
    BubblePreset
        ``(eom, pulse)`` ready for :func:`~jbubble.run_simulation`.

    References
    ----------
    Church, J. Acoust. Soc. Am. 97 (1995) 1510–1521.
    """
    eom = KellerMiksis(
        gas=PolytropicGas(gamma=gamma),
        shell=ThickShell(sigma=sigma, d_s=d_s, G_s=G_s, mu_s=mu_s),
        medium=NewtonianMedium(mu=mu),
        R0=R0,
        P_amb=P_amb,
        rho_L=rho_L,
        c_L=c_L,
    )
    pulse = ToneBurst(freq=freq, pressure=pressure, shape=Sine(), cycle_num=cycle_num)
    return BubblePreset(eom=eom, pulse=pulse)

GridSweep

Memory-efficient batched Cartesian-product parameter sweeps.

from jbubble.utils.gridsweep import GridSweep

jbubble.utils.gridsweep.GridSweep(fn, search_space, batch_size=512, progress=True)

Batched Cartesian-product sweep over named parameter axes.

Parameters:

Name Type Description Default
fn Callable

fn(**params) → PyTree Called with scalar keyword args (one per axis name); must be JAX-vmappable and jit-able.

required
search_space dict[str, ndarray]

{param_name: 1-D array of values} The Cartesian product of all axes is swept.

required
batch_size int

Number of grid points evaluated per vmapped call.

512
progress bool

Show a tqdm progress bar during iteration.

True
Source code in jbubble/utils/gridsweep.py
def __init__(
    self,
    fn: Callable,
    search_space: dict[str, jnp.ndarray],
    batch_size: int = 512,
    progress: bool = True,
) -> None:
    self.fn = fn
    self.search_space = search_space
    self.batch_size = batch_size
    self.progress = progress

    # Stable ordering so multi-index resolution is deterministic
    self._keys = sorted(search_space)
    self._axes = [jnp.asarray(search_space[k]) for k in self._keys]
    self._sizes = [int(len(a)) for a in self._axes]
    self._N = math.prod(self._sizes)

    self._eval_batch = jax.jit(jax.vmap(lambda p: self.fn(**p)))

grid_shape property

Shape of the full parameter grid (one int per axis).

total_points property

Total number of grid points.

axes property

Parameter axes in sweep order.

batches()

Lazy iterator over the grid.

Yields:

Name Type Description
params dict[str, ndarray]

Batch of parameter vectors, one per grid point in the batch.

outputs PyTree

Corresponding vmapped outputs from fn.

Source code in jbubble/utils/gridsweep.py
def batches(self):
    """Lazy iterator over the grid.

    Yields
    ------
    params : dict[str, jnp.ndarray]
        Batch of parameter vectors, one per grid point in the batch.
    outputs : PyTree
        Corresponding vmapped outputs from ``fn``.
    """
    n_batches = math.ceil(self._N / self.batch_size)
    bar = (
        tqdm(range(n_batches), desc="Grid sweep")
        if self.progress
        else range(n_batches)
    )
    for b in bar:
        start = b * self.batch_size
        end = min(start + self.batch_size, self._N)
        flat = jnp.arange(start, end)
        multi = jnp.unravel_index(flat, self._sizes)
        params = {k: self._axes[i][multi[i]] for i, k in enumerate(self._keys)}
        outputs = self._eval_batch(params)
        yield params, outputs

run()

Run the full sweep and return a grid-shaped PyTree.

Returns:

Type Description
PyTree

Each leaf has shape (*grid_shape, *leaf_shape), where grid_shape corresponds to the axes in search_space (sorted alphabetically, row-major — last axis varies fastest).

Source code in jbubble/utils/gridsweep.py
def run(self) -> PyTree:
    """Run the full sweep and return a grid-shaped PyTree.

    Returns
    -------
    PyTree
        Each leaf has shape ``(*grid_shape, *leaf_shape)``, where
        ``grid_shape`` corresponds to the axes in ``search_space``
        (sorted alphabetically, row-major — last axis varies fastest).
    """
    chunks = [outputs for _, outputs in self.batches()]
    flat = jax.tree.map(lambda *xs: jnp.concatenate(xs, axis=0), *chunks)
    shape = self.grid_shape
    return jax.tree.map(lambda x: x.reshape(*shape, *x.shape[1:]), flat)

Example

import jax.numpy as jnp
from jbubble.utils.gridsweep import GridSweep
from jbubble import run_simulation, SaveSpec
from jbubble.utils.presets import free_bubble

def peak_ratio(R0, pressure):
    preset = free_bubble(R0=R0, pressure=pressure)
    result = run_simulation(
        preset.eom, preset.pulse,
        save_spec=SaveSpec(500),
        t_max=10e-6,
    )
    return result.radius.max() / R0

sweep = GridSweep(
    fn=peak_ratio,
    search_space={
        "R0":       jnp.linspace(1e-6, 5e-6, 20),
        "pressure": jnp.array([50e3, 100e3, 200e3, 400e3]),
    },
    batch_size=256,
)

grid = sweep.run()         # shape (20, 4)
print(sweep.grid_shape)    # (20, 4)
print(sweep.total_points)  # 80

HDF5 I/O

from jbubble.utils.io import export_hdf5, load_hdf5

jbubble.utils.io.export_hdf5(path, *, metadata=None, **arrays)

Save named arrays and optional metadata to an HDF5 file.

Parameters:

Name Type Description Default
path str or Path

Output .h5 file path.

required
metadata dict

JSON-serialisable metadata stored as an attribute on the root group.

None
**arrays Any

Keyword arguments are saved as datasets. Values are converted to NumPy arrays via np.asarray.

{}

Examples:

::

result = run_simulation(eom, pulse, ...)
p_em = jax.vmap(lambda r: emission(result, r))(distances)

export_hdf5(
    "training_data.h5",
    ts=result.ts,
    R=result.state.R,
    R_dot=result.state.R_dot,
    p_emission=p_em,
    distances=distances,
    metadata={"R0": 2e-6, "freq": 1e6},
)
Source code in jbubble/utils/io.py
def export_hdf5(
    path: str | Path,
    *,
    metadata: dict[str, Any] | None = None,
    **arrays: Any,
) -> None:
    """Save named arrays and optional metadata to an HDF5 file.

    Parameters
    ----------
    path : str or Path
        Output ``.h5`` file path.
    metadata : dict, optional
        JSON-serialisable metadata stored as an attribute on the root group.
    **arrays
        Keyword arguments are saved as datasets.  Values are converted to
        NumPy arrays via ``np.asarray``.

    Examples
    --------
    ::

        result = run_simulation(eom, pulse, ...)
        p_em = jax.vmap(lambda r: emission(result, r))(distances)

        export_hdf5(
            "training_data.h5",
            ts=result.ts,
            R=result.state.R,
            R_dot=result.state.R_dot,
            p_emission=p_em,
            distances=distances,
            metadata={"R0": 2e-6, "freq": 1e6},
        )
    """
    path = Path(path)

    with h5py.File(path, "w") as f:
        for name, arr in arrays.items():
            f.create_dataset(name, data=np.asarray(arr))

        if metadata is not None:
            f.attrs["metadata"] = json.dumps(metadata)

jbubble.utils.io.load_hdf5(path)

Load arrays and metadata from an HDF5 file written by :func:export_hdf5.

Parameters:

Name Type Description Default
path str or Path

Path to the .h5 file.

required

Returns:

Name Type Description
arrays dict[str, ndarray]

All datasets in the file, keyed by name.

metadata dict

The metadata dict, or {} if none was stored.

Source code in jbubble/utils/io.py
def load_hdf5(path: str | Path) -> tuple[dict[str, np.ndarray], dict[str, Any]]:
    """Load arrays and metadata from an HDF5 file written by :func:`export_hdf5`.

    Parameters
    ----------
    path : str or Path
        Path to the ``.h5`` file.

    Returns
    -------
    arrays : dict[str, np.ndarray]
        All datasets in the file, keyed by name.
    metadata : dict
        The metadata dict, or ``{}`` if none was stored.
    """
    path = Path(path)

    with h5py.File(path, "r") as f:
        arrays = {name: np.asarray(ds) for name, ds in f.items()}
        raw = f.attrs.get("metadata", "{}")
        metadata = json.loads(raw)

    return arrays, metadata

Example

from jbubble.utils.io import export_hdf5, load_hdf5
import jax.numpy as jnp

export_hdf5(
    "results.h5",
    metadata={"description": "sweep", "freq": 1e6, "n_cycles": 5},
    R0=jnp.linspace(1e-6, 5e-6, 20),
    peak_expansion=grid,
)

arrays, meta = load_hdf5("results.h5")
print(meta["description"])          # "sweep"
print(arrays["peak_expansion"].shape)  # (20, 4)