nlsq.global_optimization.sampling¶
Sampling strategies for multi-start global optimization.
This module provides various sampling methods for generating initial parameter guesses in multi-start optimization scenarios.
Sampling Algorithms for Global Optimization¶
This module provides sampling algorithms for generating starting points in multi-start global optimization. Includes Latin Hypercube Sampling (LHS), Sobol sequences, and Halton sequences.
All samplers generate samples in the unit hypercube [0, 1]^n_dims, which can then be transformed to parameter bounds using scale_samples_to_bounds().
Key Features¶
Latin Hypercube Sampling with stratification guarantees
Sobol quasi-random sequences for low-discrepancy sampling
Halton sequences using prime bases
Bounds transformation utilities
Centering around initial parameter estimates
Examples
Generate LHS samples:
>>> from nlsq.global_optimization.sampling import latin_hypercube_sample
>>> import jax
>>> samples = latin_hypercube_sample(10, 3, rng_key=jax.random.PRNGKey(42))
>>> samples.shape
(10, 3)
Scale samples to parameter bounds:
>>> from nlsq.global_optimization.sampling import scale_samples_to_bounds
>>> import jax.numpy as jnp
>>> lb = jnp.array([0.0, -10.0])
>>> ub = jnp.array([100.0, 10.0])
>>> scaled = scale_samples_to_bounds(samples[:, :2], lb, ub)
Use the sampler factory:
>>> from nlsq.global_optimization.sampling import get_sampler
>>> sampler = get_sampler('lhs')
>>> samples = sampler(20, 4)
- nlsq.global_optimization.sampling.center_samples_around_p0(samples, p0, scale_factor, lb, ub)[source]¶
Center samples around p0 with specified scale factor.
Instead of sampling uniformly in [lb, ub], this creates samples in a region centered at p0. The region width is scale_factor * (ub - lb).
- Parameters:
samples (jax.Array) – Array of shape (n_samples, n_dims) with values in [0, 1].
p0 (jax.Array) – Center point (initial parameter guess), shape (n_dims,).
scale_factor (float) – Width of exploration region as fraction of (ub - lb). 0.5 means explore within 50% of the range centered at p0.
lb (jax.Array) – Lower bounds for each dimension, shape (n_dims,).
ub (jax.Array) – Upper bounds for each dimension, shape (n_dims,).
- Returns:
Centered samples, clipped to stay within [lb, ub].
- Return type:
Notes
The exploration region is: - Half-width = scale_factor * (ub - lb) / 2 - Region = [p0 - half_width, p0 + half_width] - Clipped to [lb, ub] to ensure valid samples
Examples
>>> import jax.numpy as jnp >>> samples = jnp.array([[0.0], [0.5], [1.0]]) >>> p0 = jnp.array([5.0]) >>> lb = jnp.array([0.0]) >>> ub = jnp.array([10.0]) >>> centered = center_samples_around_p0(samples, p0, 0.5, lb, ub) >>> # With scale_factor=0.5, half-width = 0.5 * 10 / 2 = 2.5 >>> # Region = [2.5, 7.5] >>> # Sample at 0.5 maps to center (5.0)
- nlsq.global_optimization.sampling.get_sampler(sampler_type)[source]¶
Get a sampler function by type name.
Factory function that returns the appropriate sampling function based on the sampler type string.
- Parameters:
sampler_type (str) – Type of sampler. One of: ‘lhs’, ‘sobol’, ‘halton’ (case-insensitive).
- Returns:
Sampler function with signature
(n_samples, n_dims, **kwargs) -> jax.Array.- Return type:
Callable
- Raises:
ValueError – If sampler_type is not recognized.
Examples
>>> sampler = get_sampler('lhs') >>> samples = sampler(10, 3) >>> samples.shape (10, 3)
>>> sampler = get_sampler('sobol') >>> samples = sampler(16, 4)
- nlsq.global_optimization.sampling.halton_sample(n_samples, n_dims, skip=0)[source]¶
Generate Halton quasi-random samples in the unit hypercube.
Halton sequences use different prime number bases for each dimension to create low-discrepancy samples.
- Parameters:
- Returns:
Array of shape (n_samples, n_dims) with values in [0, 1).
- Return type:
Notes
Each dimension uses a different prime base: 2, 3, 5, 7, 11, … The Halton sequence in base b for index n is computed by reflecting the base-b representation of n about the decimal point.
Examples
>>> samples = halton_sample(20, 5) >>> samples.shape (20, 5)
- nlsq.global_optimization.sampling.latin_hypercube_sample(n_samples, n_dims, rng_key=None)[source]¶
Generate Latin Hypercube samples in the unit hypercube.
Latin Hypercube Sampling divides each dimension into n_samples equal strata and ensures exactly one sample falls in each stratum for each dimension. This provides better coverage than random sampling.
- Parameters:
- Returns:
Array of shape (n_samples, n_dims) with values in [0, 1).
- Return type:
Notes
The stratification property ensures that when samples are sorted in any dimension, sample i falls in the stratum [i/n, (i+1)/n) for that dimension.
Examples
>>> import jax >>> samples = latin_hypercube_sample(10, 3, rng_key=jax.random.PRNGKey(42)) >>> samples.shape (10, 3)
>>> # Verify stratification >>> import numpy as np >>> sorted_dim0 = np.sort(np.asarray(samples[:, 0])) >>> for i, s in enumerate(sorted_dim0): ... assert i / 10 <= s < (i + 1) / 10
- nlsq.global_optimization.sampling.scale_samples_to_bounds(samples, lb, ub)[source]¶
Transform samples from [0, 1] to parameter bounds [lb, ub].
- Parameters:
- Returns:
Scaled samples with values in [lb, ub].
- Return type:
Examples
>>> import jax.numpy as jnp >>> samples = jnp.array([[0.0, 0.5], [1.0, 0.0]]) >>> lb = jnp.array([0.0, -10.0]) >>> ub = jnp.array([10.0, 10.0]) >>> scaled = scale_samples_to_bounds(samples, lb, ub) >>> # First sample: [0, 0] in [0,1] -> [0, 0] in scaled space >>> # Second sample: [1, 0] in [0,1] -> [10, -10] in scaled space
- nlsq.global_optimization.sampling.sobol_sample(n_samples, n_dims, skip=0)[source]¶
Generate Sobol quasi-random samples in the unit hypercube.
Sobol sequences are deterministic low-discrepancy sequences that provide excellent coverage of the parameter space.
- Parameters:
- Returns:
Array of shape (n_samples, n_dims) with values in [0, 1].
- Return type:
Notes
This implementation uses a Gray code approach for efficient generation. For high dimensions (>21), consider using scipy.stats.qmc.Sobol.
Examples
>>> samples = sobol_sample(16, 4) >>> samples.shape (16, 4)
Sampling Methods¶
The module supports the following sampling strategies:
Method |
Description |
|---|---|
|
Latin Hypercube Sampling - space-filling with good coverage |
|
Sobol quasi-random sequences - low discrepancy sampling |
|
Halton sequences - deterministic quasi-random sampling |
|
Uniform random sampling |
Usage Example¶
from nlsq.global_optimization.sampling import create_sampler
# Create a Latin Hypercube sampler
sampler = create_sampler("lhs", n_dims=3, bounds=[(0, 10), (0, 5), (0, 1)])
# Generate 20 starting points
points = sampler.sample(n_points=20)