Skip to content

Generators

The generators module provides three levels of synthetic InSAR data creation: single interferograms, time series, and batch generation for ML training.

Single Interferogram

generate_synthetic_insar

generate_synthetic_insar(Mw: Optional[float] = None, M0: Optional[float] = None, strike_deg: float = 0.0, dip_deg: float = 45.0, rake_deg: float = 90.0, xcen_km: float = 0.0, ycen_km: float = 0.0, depth_km: float = 10.0, grid_size: Optional[int] = None, grid_extent_km: float = 50.0, grid_spacing_km: Optional[float] = None, nu: float = DEFAULT_POISSON_RATIO, mu: float = DEFAULT_SHEAR_MODULUS_PA, satellite: Optional[str] = None, orbit: str = 'ascending', incidence_deg: Optional[float] = None, heading_deg: Optional[float] = None, wavelength_m: Optional[float] = None, add_noise: bool = True, noise_amplitude_m: float = 0.005, add_orbital_ramp: bool = False, wrap: bool = True, seed: Optional[int] = None) -> Dict

Generate synthetic InSAR data for earthquake deformation.

This is the main generator function that creates a complete synthetic interferogram including displacement field, phase, and noise. Uses the Davis (1986) point source model.

PARAMETER DESCRIPTION
Mw

Moment magnitude (provide either Mw or M0)

TYPE: float DEFAULT: None

M0

Scalar seismic moment in N·m

TYPE: float DEFAULT: None

strike_deg

Fault strike in degrees (0-360, clockwise from North)

TYPE: float DEFAULT: 0.0

dip_deg

Fault dip in degrees (0-90, from horizontal)

TYPE: float DEFAULT: 45.0

rake_deg

Slip rake in degrees (-180 to 180) - 0°: left-lateral strike-slip - 90°: thrust/reverse - ±180°: right-lateral strike-slip - -90°: normal fault

TYPE: float DEFAULT: 90.0

xcen_km

Epicenter location in km

TYPE: float DEFAULT: 0.0

ycen_km

Epicenter location in km

TYPE: float DEFAULT: 0.0

depth_km

Source depth in km (positive downward)

TYPE: float DEFAULT: 10.0

grid_size

Number of pixels for the grid (e.g., 128 for 128x128). If provided, grid_spacing_km is calculated automatically.

TYPE: int DEFAULT: None

grid_extent_km

Half-width of the grid in km (default: 50 km, so 100 km total)

TYPE: float DEFAULT: 50.0

grid_spacing_km

Grid spacing in km. If not provided and grid_size is given, calculated as: 2 * grid_extent_km / (grid_size - 1)

TYPE: float DEFAULT: None

nu

Poisson's ratio (default: 0.25)

TYPE: float DEFAULT: DEFAULT_POISSON_RATIO

mu

Shear modulus in Pa (default: 30 GPa)

TYPE: float DEFAULT: DEFAULT_SHEAR_MODULUS_PA

satellite

Satellite name for automatic geometry. Options: 'sentinel1', 'alos2', 'terrasar', 'cosmo', 'radarsat2', 'nisar', etc. If provided, overrides incidence_deg, heading_deg, wavelength_m

TYPE: str DEFAULT: None

orbit

Orbit direction: 'ascending' or 'descending' (used with satellite)

TYPE: str DEFAULT: 'ascending'

incidence_deg

Radar incidence angle in degrees (manual specification)

TYPE: float DEFAULT: None

heading_deg

Satellite heading in degrees (manual specification)

TYPE: float DEFAULT: None

wavelength_m

Radar wavelength in meters (manual specification)

TYPE: float DEFAULT: None

add_noise

Whether to add random noise

TYPE: bool DEFAULT: True

noise_amplitude_m

Noise standard deviation in meters (default: 5 mm)

TYPE: float DEFAULT: 0.005

add_orbital_ramp

Whether to add orbital ramp artifact

TYPE: bool DEFAULT: False

wrap

Whether to wrap phase to [-π, π]

TYPE: bool DEFAULT: True

seed

Random seed for reproducibility

TYPE: int DEFAULT: None

RETURNS DESCRIPTION
result

Dictionary containing: - X_km, Y_km: coordinate grids (km) - Ue, Un, Uz: displacement components (m) - los_displacement: LOS displacement (m) - phase_unwrapped: unwrapped phase (rad) - phase_wrapped: wrapped phase (rad) - phase_noisy: wrapped phase with noise (rad) - metadata: source and processing parameters

TYPE: dict

Examples:

>>> result = generate_synthetic_insar(
...     Mw=6.0,
...     strike_deg=45,
...     dip_deg=30,
...     rake_deg=90,
...     depth_km=10,
...     satellite='sentinel1',
...     orbit='ascending'
... )
Source code in src/eq_insar/generators/single.py
def generate_synthetic_insar(
    # Source parameters
    Mw: Optional[float] = None,
    M0: Optional[float] = None,
    strike_deg: float = 0.0,
    dip_deg: float = 45.0,
    rake_deg: float = 90.0,
    # Source location
    xcen_km: float = 0.0,
    ycen_km: float = 0.0,
    depth_km: float = 10.0,
    # Grid parameters
    grid_size: Optional[int] = None,
    grid_extent_km: float = 50.0,
    grid_spacing_km: Optional[float] = None,
    # Elastic parameters
    nu: float = DEFAULT_POISSON_RATIO,
    mu: float = DEFAULT_SHEAR_MODULUS_PA,
    # InSAR geometry - can specify satellite or manual
    satellite: Optional[str] = None,
    orbit: str = "ascending",
    incidence_deg: Optional[float] = None,
    heading_deg: Optional[float] = None,
    wavelength_m: Optional[float] = None,
    # Noise parameters
    add_noise: bool = True,
    noise_amplitude_m: float = 0.005,
    add_orbital_ramp: bool = False,
    # Output options
    wrap: bool = True,
    seed: Optional[int] = None,
) -> Dict:
    """
    Generate synthetic InSAR data for earthquake deformation.

    This is the main generator function that creates a complete synthetic
    interferogram including displacement field, phase, and noise.
    Uses the Davis (1986) point source model.

    Parameters
    ----------
    Mw : float, optional
        Moment magnitude (provide either Mw or M0)
    M0 : float, optional
        Scalar seismic moment in N·m
    strike_deg : float
        Fault strike in degrees (0-360, clockwise from North)
    dip_deg : float
        Fault dip in degrees (0-90, from horizontal)
    rake_deg : float
        Slip rake in degrees (-180 to 180)
        - 0°: left-lateral strike-slip
        - 90°: thrust/reverse
        - ±180°: right-lateral strike-slip
        - -90°: normal fault
    xcen_km, ycen_km : float
        Epicenter location in km
    depth_km : float
        Source depth in km (positive downward)

    grid_size : int, optional
        Number of pixels for the grid (e.g., 128 for 128x128).
        If provided, grid_spacing_km is calculated automatically.
    grid_extent_km : float
        Half-width of the grid in km (default: 50 km, so 100 km total)
    grid_spacing_km : float, optional
        Grid spacing in km. If not provided and grid_size is given,
        calculated as: 2 * grid_extent_km / (grid_size - 1)

    nu : float
        Poisson's ratio (default: 0.25)
    mu : float
        Shear modulus in Pa (default: 30 GPa)

    satellite : str, optional
        Satellite name for automatic geometry. Options:
        'sentinel1', 'alos2', 'terrasar', 'cosmo', 'radarsat2', 'nisar', etc.
        If provided, overrides incidence_deg, heading_deg, wavelength_m
    orbit : str
        Orbit direction: 'ascending' or 'descending' (used with satellite)
    incidence_deg : float, optional
        Radar incidence angle in degrees (manual specification)
    heading_deg : float, optional
        Satellite heading in degrees (manual specification)
    wavelength_m : float, optional
        Radar wavelength in meters (manual specification)

    add_noise : bool
        Whether to add random noise
    noise_amplitude_m : float
        Noise standard deviation in meters (default: 5 mm)
    add_orbital_ramp : bool
        Whether to add orbital ramp artifact

    wrap : bool
        Whether to wrap phase to [-π, π]
    seed : int, optional
        Random seed for reproducibility

    Returns
    -------
    result : dict
        Dictionary containing:
        - X_km, Y_km: coordinate grids (km)
        - Ue, Un, Uz: displacement components (m)
        - los_displacement: LOS displacement (m)
        - phase_unwrapped: unwrapped phase (rad)
        - phase_wrapped: wrapped phase (rad)
        - phase_noisy: wrapped phase with noise (rad)
        - metadata: source and processing parameters

    Examples
    --------
    >>> result = generate_synthetic_insar(
    ...     Mw=6.0,
    ...     strike_deg=45,
    ...     dip_deg=30,
    ...     rake_deg=90,
    ...     depth_km=10,
    ...     satellite='sentinel1',
    ...     orbit='ascending'
    ... )
    """
    # Validate source parameters
    _validate_source_parameters(Mw, M0, strike_deg, dip_deg, rake_deg, depth_km)

    if seed is not None:
        np.random.seed(seed)

    # Handle grid_size parameter
    if grid_size is not None:
        # Calculate grid_spacing_km from grid_size and grid_extent_km
        grid_spacing_km = (2 * grid_extent_km) / (grid_size - 1)
    elif grid_spacing_km is None:
        # Default spacing if neither grid_size nor grid_spacing_km provided
        grid_spacing_km = 0.5

    # Get seismic moment
    if M0 is None:
        M0 = mw_to_m0(Mw)
    else:
        Mw = m0_to_mw(M0)

    # Get satellite configuration or use manual parameters
    if satellite is not None:
        sat = get_satellite(satellite)
        if incidence_deg is None:
            incidence_deg = sat.incidence_deg
        if heading_deg is None:
            heading_deg = sat.get_heading(orbit)
        if wavelength_m is None:
            wavelength_m = sat.wavelength_m
        satellite_name = sat.name
    else:
        # Defaults if nothing specified
        if incidence_deg is None:
            incidence_deg = 33.0
        if heading_deg is None:
            heading_deg = -13.0
        if wavelength_m is None:
            wavelength_m = 0.05546  # Sentinel-1
        satellite_name = "custom"

    # Create coordinate grids
    if grid_size is not None:
        # Use linspace for exact grid_size
        xs_km = np.linspace(-grid_extent_km, grid_extent_km, grid_size)
        ys_km = np.linspace(-grid_extent_km, grid_extent_km, grid_size)
    else:
        # Use arange with spacing
        xs_km = np.arange(-grid_extent_km, grid_extent_km + grid_spacing_km, grid_spacing_km)
        ys_km = np.arange(-grid_extent_km, grid_extent_km + grid_spacing_km, grid_spacing_km)
    X_km, Y_km = np.meshgrid(xs_km, ys_km)

    # Convert to meters for physics calculations
    X_m = X_km * 1000.0
    Y_m = Y_km * 1000.0
    xcen_m = xcen_km * 1000.0
    ycen_m = ycen_km * 1000.0
    depth_m = depth_km * 1000.0

    # Build moment tensor
    Mxx, Myy, Mzz, Mxy, Myz, Mzx = double_couple_moment_tensor(
        strike_deg, dip_deg, rake_deg, M0
    )

    # Compute displacement using Davis point source
    Ue, Un, Uz = davis_point_source(
        X_m, Y_m, xcen_m, ycen_m, depth_m,
        Mxx, Myy, Mzz, Mxy, Myz, Mzx,
        nu=nu, mu=mu
    )

    # Compute LOS displacement
    d_los = compute_los_displacement(
        Ue, Un, Uz,
        incidence_deg=incidence_deg,
        heading_deg=heading_deg
    )

    # Convert to phase
    phase_unwrapped = displacement_to_phase(d_los, wavelength_m)
    phase_wrapped = wrap_phase(phase_unwrapped) if wrap else phase_unwrapped

    # Generate noise
    ny, nx = X_km.shape

    # Initialize noisy phase
    phase_noisy = phase_unwrapped.copy()

    if add_noise:
        # Random noise
        noise = generate_random_noise(
            (ny, nx),
            amplitude_m=noise_amplitude_m,
            seed=seed
        )
        noise_phase = displacement_to_phase(noise, wavelength_m)
        phase_noisy = phase_noisy + noise_phase

        # Orbital ramp
        if add_orbital_ramp:
            ramp = generate_orbital_ramp(
                (ny, nx),
                pixel_size_km=grid_spacing_km,
                seed=seed + 300 if seed else None
            )
            ramp_phase = displacement_to_phase(ramp, wavelength_m)
            phase_noisy = phase_noisy + ramp_phase

    # Wrap noisy phase
    if wrap:
        phase_noisy = wrap_phase(phase_noisy)

    # Package metadata
    metadata = {
        # Source parameters
        "Mw": Mw,
        "M0_Nm": M0,
        "strike_deg": strike_deg,
        "dip_deg": dip_deg,
        "rake_deg": rake_deg,
        "xcen_km": xcen_km,
        "ycen_km": ycen_km,
        "depth_km": depth_km,
        "source_type": "davis",
        # Grid parameters
        "grid_extent_km": grid_extent_km,
        "grid_spacing_km": grid_spacing_km,
        # Elastic parameters
        "nu": nu,
        "mu_Pa": mu,
        # InSAR geometry
        "satellite": satellite_name,
        "orbit": orbit,
        "incidence_deg": incidence_deg,
        "heading_deg": heading_deg,
        "wavelength_m": wavelength_m,
        # Noise parameters
        "noise_amplitude_m": noise_amplitude_m,
    }

    return {
        "X_km": X_km,
        "Y_km": Y_km,
        "Ue": Ue,
        "Un": Un,
        "Uz": Uz,
        "los_displacement": d_los,
        "phase_unwrapped": phase_unwrapped,
        "phase_wrapped": phase_wrapped,
        "phase_noisy": phase_noisy,
        "metadata": metadata,
    }

Time Series

generate_timeseries

generate_timeseries(Mw: Optional[float] = None, M0: Optional[float] = None, strike_deg: float = 0.0, dip_deg: float = 45.0, rake_deg: float = 90.0, xcen_km: float = 0.0, ycen_km: float = 0.0, depth_km: float = 10.0, grid_size: Optional[int] = None, grid_extent_km: float = 50.0, grid_spacing_km: Optional[float] = None, nu: float = DEFAULT_POISSON_RATIO, mu: float = DEFAULT_SHEAR_MODULUS_PA, satellite: Optional[str] = None, orbit: str = 'ascending', incidence_deg: Optional[float] = None, heading_deg: Optional[float] = None, wavelength_m: Optional[float] = None, n_pre: int = 5, n_event: int = 1, n_post: int = 5, noise_amplitude_m: float = 0.005, wrap: bool = True, output_type: str = 'phase', deformation_threshold_m: float = 0.005, seed: Optional[int] = None) -> Dict

Generate time series of synthetic InSAR data.

Creates a sequence with: - Pre-event frames: noise only (no deformation signal) - Event frames: deformation signal + noise - Post-event frames: noise only

This format is designed for training ML models to detect and segment earthquake deformation signals.

PARAMETER DESCRIPTION
Mw

Moment magnitude

TYPE: float DEFAULT: None

M0

Seismic moment in N·m

TYPE: float DEFAULT: None

strike_deg

Fault geometry parameters

TYPE: float DEFAULT: 0.0

dip_deg

Fault geometry parameters

TYPE: float DEFAULT: 0.0

rake_deg

Fault geometry parameters

TYPE: float DEFAULT: 0.0

xcen_km

Epicenter location in km

TYPE: float DEFAULT: 0.0

ycen_km

Epicenter location in km

TYPE: float DEFAULT: 0.0

depth_km

Source depth in km

TYPE: float DEFAULT: 10.0

grid_size

Number of pixels for the grid (e.g., 128 for 128x128). If provided, grid_spacing_km is calculated automatically.

TYPE: int DEFAULT: None

grid_extent_km

Half-width of the grid in km (default: 50 km)

TYPE: float DEFAULT: 50.0

grid_spacing_km

Grid spacing in km. If not provided and grid_size is given, calculated as: 2 * grid_extent_km / (grid_size - 1)

TYPE: float DEFAULT: None

satellite

Satellite configuration ('sentinel1', 'alos2', etc.)

TYPE: str DEFAULT: None

orbit

'ascending' or 'descending'

TYPE: str DEFAULT: 'ascending'

n_pre

Number of pre-event frames (noise only)

TYPE: int DEFAULT: 5

n_event

Number of event frames (signal + noise)

TYPE: int DEFAULT: 1

n_post

Number of post-event frames (noise only)

TYPE: int DEFAULT: 5

noise_amplitude_m

Noise standard deviation in meters

TYPE: float DEFAULT: 0.005

wrap

Wrap phase to [-π, π]

TYPE: bool DEFAULT: True

output_type

'phase' (default), 'displacement', or 'both'

TYPE: str DEFAULT: 'phase'

deformation_threshold_m

Threshold for creating binary deformation labels

TYPE: float DEFAULT: 0.005

seed

Random seed for reproducibility

TYPE: int DEFAULT: None

RETURNS DESCRIPTION
result

Dictionary containing: - All fields from generate_synthetic_insar (static deformation) - timeseries: (n_total, ny, nx) array of frames - labels: (n_total, ny, nx) binary segmentation masks - metadata: includes n_pre, n_event, n_post

TYPE: dict

Examples:

>>> result = generate_timeseries(
...     Mw=6.0,
...     satellite='sentinel1',
...     n_pre=5,
...     n_event=1,
...     n_post=5
... )
>>> X = result['timeseries']  # Shape: (11, ny, nx)
>>> y = result['labels']      # Binary masks
Source code in src/eq_insar/generators/timeseries.py
def generate_timeseries(
    # Source parameters
    Mw: Optional[float] = None,
    M0: Optional[float] = None,
    strike_deg: float = 0.0,
    dip_deg: float = 45.0,
    rake_deg: float = 90.0,
    xcen_km: float = 0.0,
    ycen_km: float = 0.0,
    depth_km: float = 10.0,
    # Grid parameters
    grid_size: Optional[int] = None,
    grid_extent_km: float = 50.0,
    grid_spacing_km: Optional[float] = None,
    # Elastic parameters
    nu: float = DEFAULT_POISSON_RATIO,
    mu: float = DEFAULT_SHEAR_MODULUS_PA,
    # InSAR geometry
    satellite: Optional[str] = None,
    orbit: str = "ascending",
    incidence_deg: Optional[float] = None,
    heading_deg: Optional[float] = None,
    wavelength_m: Optional[float] = None,
    # Time series parameters
    n_pre: int = 5,
    n_event: int = 1,
    n_post: int = 5,
    # Noise parameters
    noise_amplitude_m: float = 0.005,
    # Output options
    wrap: bool = True,
    output_type: str = "phase",  # 'phase', 'displacement', 'both'
    deformation_threshold_m: float = 0.005,  # 5mm for label mask
    seed: Optional[int] = None,
) -> Dict:
    """
    Generate time series of synthetic InSAR data.

    Creates a sequence with:
    - Pre-event frames: noise only (no deformation signal)
    - Event frames: deformation signal + noise
    - Post-event frames: noise only

    This format is designed for training ML models to detect and
    segment earthquake deformation signals.

    Parameters
    ----------
    Mw : float, optional
        Moment magnitude
    M0 : float, optional
        Seismic moment in N·m
    strike_deg, dip_deg, rake_deg : float
        Fault geometry parameters
    xcen_km, ycen_km : float
        Epicenter location in km
    depth_km : float
        Source depth in km

    grid_size : int, optional
        Number of pixels for the grid (e.g., 128 for 128x128).
        If provided, grid_spacing_km is calculated automatically.
    grid_extent_km : float
        Half-width of the grid in km (default: 50 km)
    grid_spacing_km : float, optional
        Grid spacing in km. If not provided and grid_size is given,
        calculated as: 2 * grid_extent_km / (grid_size - 1)

    satellite : str, optional
        Satellite configuration ('sentinel1', 'alos2', etc.)
    orbit : str
        'ascending' or 'descending'

    n_pre : int
        Number of pre-event frames (noise only)
    n_event : int
        Number of event frames (signal + noise)
    n_post : int
        Number of post-event frames (noise only)

    noise_amplitude_m : float
        Noise standard deviation in meters

    wrap : bool
        Wrap phase to [-π, π]
    output_type : str
        'phase' (default), 'displacement', or 'both'
    deformation_threshold_m : float
        Threshold for creating binary deformation labels
    seed : int, optional
        Random seed for reproducibility

    Returns
    -------
    result : dict
        Dictionary containing:
        - All fields from generate_synthetic_insar (static deformation)
        - timeseries: (n_total, ny, nx) array of frames
        - labels: (n_total, ny, nx) binary segmentation masks
        - metadata: includes n_pre, n_event, n_post

    Examples
    --------
    >>> result = generate_timeseries(
    ...     Mw=6.0,
    ...     satellite='sentinel1',
    ...     n_pre=5,
    ...     n_event=1,
    ...     n_post=5
    ... )
    >>> X = result['timeseries']  # Shape: (11, ny, nx)
    >>> y = result['labels']      # Binary masks
    """
    if seed is not None:
        np.random.seed(seed)

    # Generate the static deformation signal (without noise)
    result = generate_synthetic_insar(
        Mw=Mw,
        M0=M0,
        strike_deg=strike_deg,
        dip_deg=dip_deg,
        rake_deg=rake_deg,
        xcen_km=xcen_km,
        ycen_km=ycen_km,
        depth_km=depth_km,
        grid_size=grid_size,
        grid_extent_km=grid_extent_km,
        grid_spacing_km=grid_spacing_km,
        nu=nu,
        mu=mu,
        satellite=satellite,
        orbit=orbit,
        incidence_deg=incidence_deg,
        heading_deg=heading_deg,
        wavelength_m=wavelength_m,
        add_noise=False,  # We'll add noise per-frame
        wrap=False,
        seed=seed,
    )

    ny, nx = result["X_km"].shape
    n_total = n_pre + n_event + n_post

    # Get wavelength for phase conversion
    wl = result["metadata"]["wavelength_m"]

    # Initialize output arrays
    timeseries = np.zeros((n_total, ny, nx), dtype=np.float32)

    # Create deformation labels (binary mask)
    labels = np.zeros((n_total, ny, nx), dtype=np.int32)
    deformation_mask = np.abs(result["los_displacement"]) > deformation_threshold_m

    # Generate frames
    for t in range(n_total):
        # Generate per-frame random noise (different each frame)
        noise = generate_random_noise(
            (ny, nx),
            amplitude_m=noise_amplitude_m,
            seed=seed + t * 100 if seed else None
        )

        # Check if this is an event frame
        is_event_frame = n_pre <= t < n_pre + n_event

        if is_event_frame:
            # Event frame: signal + noise
            if output_type == "displacement":
                timeseries[t] = result["los_displacement"] + noise
            else:  # phase
                signal_phase = result["phase_unwrapped"]
                noise_phase = displacement_to_phase(noise, wl)
                total_phase = signal_phase + noise_phase
                timeseries[t] = wrap_phase(total_phase) if wrap else total_phase

            # Set label
            labels[t] = deformation_mask.astype(np.int32)
        else:
            # Non-event frame: noise only
            if output_type == "displacement":
                timeseries[t] = noise
            else:  # phase
                noise_phase = displacement_to_phase(noise, wl)
                timeseries[t] = wrap_phase(noise_phase) if wrap else noise_phase
            # Labels remain 0

    # Add time series data to result
    result["timeseries"] = timeseries
    result["labels"] = labels

    # Update metadata
    result["metadata"]["n_pre"] = n_pre
    result["metadata"]["n_event"] = n_event
    result["metadata"]["n_post"] = n_post
    result["metadata"]["n_total"] = n_total
    result["metadata"]["output_type"] = output_type
    result["metadata"]["deformation_threshold_m"] = deformation_threshold_m

    return result

Batch Generation

generate_training_batch

generate_training_batch(n_samples: int, grid_size: Optional[int] = None, grid_extent_km: float = 50.0, grid_spacing_km: Optional[float] = None, mw_range: Tuple[float, float] = (4.5, 7.0), depth_range_km: Tuple[float, float] = (5.0, 20.0), fault_type: Optional[str] = None, satellite: str = 'sentinel1', orbit: str = 'ascending', n_pre: int = 3, n_event: int = 1, n_post: int = 3, noise_range_m: Tuple[float, float] = (0.002, 0.008), wrap: bool = True, output_type: str = 'phase', seed: Optional[int] = None, verbose: bool = True) -> List[Dict]

Generate a batch of training samples with random parameters.

Creates multiple synthetic InSAR time series with varied earthquake parameters, suitable for training ML models.

PARAMETER DESCRIPTION
n_samples

Number of samples to generate

TYPE: int

grid_size

Number of pixels for the grid (e.g., 128 for 128x128). If provided, grid_spacing_km is calculated automatically.

TYPE: int DEFAULT: None

grid_extent_km

Half-width of the grid in km (default: 50 km)

TYPE: float DEFAULT: 50.0

grid_spacing_km

Grid spacing in km. If not provided and grid_size is given, calculated as: 2 * grid_extent_km / (grid_size - 1)

TYPE: float DEFAULT: None

mw_range

Range of moment magnitudes

TYPE: tuple(min, max) DEFAULT: (4.5, 7.0)

depth_range_km

Range of depths in km

TYPE: tuple(min, max) DEFAULT: (5.0, 20.0)

fault_type

Constrain fault type: 'strike-slip', 'thrust', 'normal', or None

TYPE: str DEFAULT: None

satellite

Satellite configuration name

TYPE: str DEFAULT: 'sentinel1'

orbit

'ascending' or 'descending'

TYPE: str DEFAULT: 'ascending'

n_pre

Time series structure

TYPE: int DEFAULT: 3

n_event

Time series structure

TYPE: int DEFAULT: 3

n_post

Time series structure

TYPE: int DEFAULT: 3

noise_range_m

Range of noise amplitudes

TYPE: tuple(min, max) DEFAULT: (0.002, 0.008)

wrap

Whether to wrap phase

TYPE: bool DEFAULT: True

output_type

'phase' or 'displacement'

TYPE: str DEFAULT: 'phase'

seed

Random seed for reproducibility

TYPE: int DEFAULT: None

verbose

Print progress messages

TYPE: bool DEFAULT: True

RETURNS DESCRIPTION
samples

List of sample dictionaries from generate_timeseries

TYPE: list of dict

Examples:

>>> batch = generate_training_batch(
...     n_samples=100,
...     mw_range=(5.0, 6.5),
...     satellite='sentinel1'
... )
>>> X = np.stack([s['timeseries'] for s in batch])  # (100, T, H, W)
>>> y = np.stack([s['labels'] for s in batch])      # (100, T, H, W)
Source code in src/eq_insar/generators/batch.py
def generate_training_batch(
    n_samples: int,
    # Grid parameters
    grid_size: Optional[int] = None,
    grid_extent_km: float = 50.0,
    grid_spacing_km: Optional[float] = None,
    # Source parameters
    mw_range: Tuple[float, float] = (4.5, 7.0),
    depth_range_km: Tuple[float, float] = (5.0, 20.0),
    fault_type: Optional[str] = None,
    # InSAR geometry
    satellite: str = "sentinel1",
    orbit: str = "ascending",
    # Time series parameters
    n_pre: int = 3,
    n_event: int = 1,
    n_post: int = 3,
    # Noise parameters
    noise_range_m: Tuple[float, float] = (0.002, 0.008),
    # Output options
    wrap: bool = True,
    output_type: str = "phase",
    seed: Optional[int] = None,
    verbose: bool = True,
) -> List[Dict]:
    """
    Generate a batch of training samples with random parameters.

    Creates multiple synthetic InSAR time series with varied earthquake
    parameters, suitable for training ML models.

    Parameters
    ----------
    n_samples : int
        Number of samples to generate
    grid_size : int, optional
        Number of pixels for the grid (e.g., 128 for 128x128).
        If provided, grid_spacing_km is calculated automatically.
    grid_extent_km : float
        Half-width of the grid in km (default: 50 km)
    grid_spacing_km : float, optional
        Grid spacing in km. If not provided and grid_size is given,
        calculated as: 2 * grid_extent_km / (grid_size - 1)
    mw_range : tuple (min, max)
        Range of moment magnitudes
    depth_range_km : tuple (min, max)
        Range of depths in km
    fault_type : str, optional
        Constrain fault type: 'strike-slip', 'thrust', 'normal', or None
    satellite : str
        Satellite configuration name
    orbit : str
        'ascending' or 'descending'
    n_pre, n_event, n_post : int
        Time series structure
    noise_range_m : tuple (min, max)
        Range of noise amplitudes
    wrap : bool
        Whether to wrap phase
    output_type : str
        'phase' or 'displacement'
    seed : int, optional
        Random seed for reproducibility
    verbose : bool
        Print progress messages

    Returns
    -------
    samples : list of dict
        List of sample dictionaries from generate_timeseries

    Examples
    --------
    >>> batch = generate_training_batch(
    ...     n_samples=100,
    ...     mw_range=(5.0, 6.5),
    ...     satellite='sentinel1'
    ... )
    >>> X = np.stack([s['timeseries'] for s in batch])  # (100, T, H, W)
    >>> y = np.stack([s['labels'] for s in batch])      # (100, T, H, W)
    """
    samples = []

    for i in range(n_samples):
        sample_seed = seed + i if seed else None

        # Sample earthquake parameters
        params = sample_earthquake_parameters(
            mw_range=mw_range,
            depth_range_km=depth_range_km,
            location_range_km=grid_extent_km * 0.6,
            fault_type=fault_type,
            seed=sample_seed,
        )

        # Sample noise amplitude
        if sample_seed is not None:
            np.random.seed(sample_seed + 1000)
        noise_amplitude = np.random.uniform(*noise_range_m)

        # Generate time series
        result = generate_timeseries(
            **params,
            grid_size=grid_size,
            grid_extent_km=grid_extent_km,
            grid_spacing_km=grid_spacing_km,
            satellite=satellite,
            orbit=orbit,
            n_pre=n_pre,
            n_event=n_event,
            n_post=n_post,
            noise_amplitude_m=noise_amplitude,
            wrap=wrap,
            output_type=output_type,
            seed=sample_seed,
        )

        samples.append(result)

        if verbose and (i + 1) % max(1, n_samples // 10) == 0:
            print(f"Generated {i + 1}/{n_samples} samples")

    return samples

batch_to_arrays

batch_to_arrays(batch: List[Dict], include_metadata: bool = False) -> Dict

Convert batch of samples to stacked numpy arrays for ML training.

PARAMETER DESCRIPTION
batch

Output from generate_training_batch

TYPE: list of dict

include_metadata

Include metadata as a list

TYPE: bool DEFAULT: False

RETURNS DESCRIPTION
data

Dictionary containing: - X: (N, T, H, W) timeseries array - y: (N, T, H, W) labels array - los: (N, H, W) static LOS displacement - metadata: list of metadata dicts (if include_metadata=True)

TYPE: dict

Source code in src/eq_insar/generators/batch.py
def batch_to_arrays(
    batch: List[Dict],
    include_metadata: bool = False,
) -> Dict:
    """
    Convert batch of samples to stacked numpy arrays for ML training.

    Parameters
    ----------
    batch : list of dict
        Output from generate_training_batch
    include_metadata : bool
        Include metadata as a list

    Returns
    -------
    data : dict
        Dictionary containing:
        - X: (N, T, H, W) timeseries array
        - y: (N, T, H, W) labels array
        - los: (N, H, W) static LOS displacement
        - metadata: list of metadata dicts (if include_metadata=True)
    """
    X = np.stack([s["timeseries"] for s in batch])
    y = np.stack([s["labels"] for s in batch])
    los = np.stack([s["los_displacement"] for s in batch])

    result = {
        "X": X,
        "y": y,
        "los_displacement": los,
    }

    if include_metadata:
        result["metadata"] = [s["metadata"] for s in batch]

    return result

Parameter Sampling

sample_earthquake_parameters

sample_earthquake_parameters(mw_range: Tuple[float, float] = (4.5, 7.0), depth_range_km: Tuple[float, float] = (5.0, 20.0), location_range_km: float = 30.0, fault_type: Optional[str] = None, seed: Optional[int] = None) -> Dict

Sample random earthquake source parameters.

Generates physically reasonable earthquake parameters with optional constraints on fault type.

PARAMETER DESCRIPTION
mw_range

Range of moment magnitudes to sample from

TYPE: tuple(min, max) DEFAULT: (4.5, 7.0)

depth_range_km

Range of depths in km to sample from

TYPE: tuple(min, max) DEFAULT: (5.0, 20.0)

location_range_km

Maximum offset from center for epicenter location

TYPE: float DEFAULT: 30.0

fault_type

Constrain fault type: 'strike-slip', 'thrust', 'normal', or None (random)

TYPE: str DEFAULT: None

seed

Random seed for reproducibility

TYPE: int DEFAULT: None

RETURNS DESCRIPTION
params

Dictionary with keys: Mw, strike_deg, dip_deg, rake_deg, xcen_km, ycen_km, depth_km

TYPE: dict

Source code in src/eq_insar/generators/batch.py
def sample_earthquake_parameters(
    mw_range: Tuple[float, float] = (4.5, 7.0),
    depth_range_km: Tuple[float, float] = (5.0, 20.0),
    location_range_km: float = 30.0,
    fault_type: Optional[str] = None,
    seed: Optional[int] = None,
) -> Dict:
    """
    Sample random earthquake source parameters.

    Generates physically reasonable earthquake parameters with optional
    constraints on fault type.

    Parameters
    ----------
    mw_range : tuple (min, max)
        Range of moment magnitudes to sample from
    depth_range_km : tuple (min, max)
        Range of depths in km to sample from
    location_range_km : float
        Maximum offset from center for epicenter location
    fault_type : str, optional
        Constrain fault type: 'strike-slip', 'thrust', 'normal', or None (random)
    seed : int, optional
        Random seed for reproducibility

    Returns
    -------
    params : dict
        Dictionary with keys: Mw, strike_deg, dip_deg, rake_deg,
        xcen_km, ycen_km, depth_km
    """
    # Validate ranges
    if mw_range[0] > mw_range[1]:
        raise ValueError(f"mw_range min ({mw_range[0]}) > max ({mw_range[1]})")
    if depth_range_km[0] > depth_range_km[1]:
        raise ValueError(
            f"depth_range_km min ({depth_range_km[0]}) > max ({depth_range_km[1]})"
        )
    if depth_range_km[0] <= 0:
        raise ValueError(
            f"depth_range_km min ({depth_range_km[0]}) must be positive"
        )

    if seed is not None:
        np.random.seed(seed)

    # Sample magnitude
    Mw = np.random.uniform(*mw_range)

    # Sample location
    xcen_km = np.random.uniform(-location_range_km, location_range_km)
    ycen_km = np.random.uniform(-location_range_km, location_range_km)
    depth_km = np.random.uniform(*depth_range_km)

    # Sample fault geometry based on type
    strike_deg = np.random.uniform(0, 360)

    if fault_type is None:
        # Random fault type
        fault_type = np.random.choice(["strike-slip", "thrust", "normal"])

    if fault_type == "strike-slip":
        dip_deg = np.random.uniform(75, 90)  # Near-vertical
        rake_deg = np.random.choice([-1, 1]) * np.random.uniform(0, 30)  # Near-horizontal slip
    elif fault_type == "thrust" or fault_type == "reverse":
        dip_deg = np.random.uniform(15, 50)  # Shallow to moderate dip
        rake_deg = np.random.uniform(60, 120)  # Dominantly thrust
    elif fault_type == "normal":
        dip_deg = np.random.uniform(45, 70)  # Moderate to steep dip
        rake_deg = np.random.uniform(-120, -60)  # Dominantly normal
    else:
        # Fully random
        dip_deg = np.random.uniform(10, 80)
        rake_deg = np.random.uniform(-180, 180)

    return {
        "Mw": Mw,
        "strike_deg": strike_deg,
        "dip_deg": dip_deg,
        "rake_deg": rake_deg,
        "xcen_km": xcen_km,
        "ycen_km": ycen_km,
        "depth_km": depth_km,
    }