Skip to content

acoustics

The jbubble.acoustics module provides models for computing the acoustic pressure radiated by a bubble, given a solved trajectory.

Emission models are not Property subclasses (they depend on the full trajectory, not just instantaneous state), and not part of the EoM (they are pure post-processing). They take a SimulationResult and a field-point distance r and return the pressure time series.

from jbubble.acoustics import IncompressibleMonopole

emission = IncompressibleMonopole(rho_L=998.0)
p_rad = emission(result, r=1e-2)   # pressure [Pa] at 1 cm, shape (num_samples,)

For multiple field-point distances, use jax.vmap:

import jax
import jax.numpy as jnp

distances = jnp.array([1e-3, 5e-3, 1e-2, 5e-2])
p_all = jax.vmap(lambda r: emission(result, r))(distances)
# shape: (4, num_samples)

jbubble.acoustics.emission.EmissionModel

Bases: Module, ABC

Acoustic emission model: bubble trajectory → radiated pressure.

Subclasses implement __call__ which takes a solved :class:~jbubble.simulation.SimulationResult and a field-point distance r and returns the radiated pressure time series.

Multiple field-point distances are handled naturally via jax.vmap::

distances = jnp.array([0.001, 0.005, 0.01])
p_all = jax.vmap(lambda r: model(result, r))(distances)
# shape (3, N)

__call__(result, r) abstractmethod

Compute radiated pressure at distance r.

Parameters:

Name Type Description Default
result SimulationResult

Solved bubble trajectory (state, state_dot, ts).

required
r float or Array

Distance from the bubble centre to the field point [m].

required

Returns:

Type Description
(Array, shape(N))

Radiated pressure [Pa] at each saved time point.

Source code in jbubble/acoustics/emission.py
@abc.abstractmethod
def __call__(
    self,
    result: SimulationResult,
    r: ArrayLike,
) -> jax.Array:
    """Compute radiated pressure at distance *r*.

    Parameters
    ----------
    result : SimulationResult
        Solved bubble trajectory (``state``, ``state_dot``, ``ts``).
    r : float or jax.Array
        Distance from the bubble centre to the field point [m].

    Returns
    -------
    jax.Array, shape (N,)
        Radiated pressure [Pa] at each saved time point.
    """
    ...

jbubble.acoustics.emission.IncompressibleMonopole

Bases: EmissionModel

Incompressible monopole radiation.

Assumes an incompressible surrounding liquid so that the radiated pressure at distance r is given by the time derivative of the volume flux:

::

p_rad(r, t) = rho_L / r · d/dt(R² Ṙ)
             = rho_L / r · (2 R Ṙ² + R² R̈)

This is the simplest acoustic emission model and is accurate when the bubble-wall Mach number Ṁ = Ṙ / c_L ≪ 1 and the field point is in the geometric near-field (r ≪ c_L / f).

Fields

rho_L : float or jax.Array Liquid density [kg/m³].

jbubble.acoustics.emission.QuasiAcoustic

Bases: EmissionModel

Quasi-acoustic emission with retarded-time correction.

Accounts for the finite speed of sound by evaluating the bubble-wall quantities at the retarded time t_ret = t − r / c_L:

::

p_rad(r, t) = rho_L R²(t_ret) / r
              · [R̈(t_ret) + 2 Ṙ²(t_ret) / R(t_ret)]

Uses linear interpolation (jnp.interp) to evaluate the trajectory at retarded times. For field points where t_ret < t_start the values are clamped to the initial (equilibrium) state — physically reasonable since the bubble is quiescent before excitation.

Fields

rho_L : float or jax.Array Liquid density [kg/m³]. c_L : float or jax.Array Speed of sound in the liquid [m/s].


Choosing an emission model

Model Accuracy Validity Speed
IncompressibleMonopole Low \(r \gg R\), \(\text{Ma} \ll 1\) Fastest — instant post-processing
QuasiAcoustic Medium \(r \gg R\), moderate Ma Fast — uses jnp.interp

The fully compressible acoustic emission model (analogous to APECSS's wave equation integration) is planned but not yet implemented. For most research purposes, QuasiAcoustic at a field point distance \(r \gg R_0\) is adequate.