att.embedding¶
Phase-space reconstruction via time-delay embedding.
- class att.embedding.TakensEmbedder(delay='auto', dimension='auto')[source]¶
Bases:
objectReconstruct a phase-space attractor from a scalar time series.
- Parameters:
- fit(X)[source]¶
Estimate parameters from data. Stores .delay_ and .dimension_.
- Parameters:
X (ndarray)
- Return type:
- class att.embedding.JointEmbedder(delays='auto', dimensions='auto')[source]¶
Bases:
objectConstruct joint delay embeddings with per-channel delay estimation.
Using “auto” is strongly recommended for systems with different timescales.
- Parameters:
- fit(channels)[source]¶
Estimate per-channel parameters. Stores .delays_ and .dimensions_.
- Parameters:
- Return type:
- transform(channels)[source]¶
Construct joint delay vectors by concatenating per-channel embeddings.
Input: list of 1D arrays, each (n_samples,) Output: (n_valid_samples, sum(dimensions))
- att.embedding.estimate_delay(X, method='ami', max_lag=100)[source]¶
Estimate optimal time delay for Takens embedding.
Uses first minimum of Average Mutual Information (Fraser & Swinney, 1986).
- Parameters:
X ((n_samples,) array)
method ("ami" (only option currently))
max_lag (maximum lag to consider)
- Returns:
int
- Return type:
optimal delay (first AMI minimum, or max_lag if no minimum found)
- att.embedding.estimate_dimension(X, delay, method='fnn', max_dim=10, threshold=0.01, rtol=15.0, atol=2.0)[source]¶
Estimate minimal embedding dimension via False Nearest Neighbors.
- Parameters:
X ((n_samples,) 1D time series)
delay (time delay for embedding)
method ("fnn" (only option currently))
max_dim (maximum dimension to test)
threshold (FNN fraction below which to stop (default 0.01 = 1%))
rtol (relative distance threshold for FNN criterion 1)
atol (absolute distance threshold for FNN criterion 2)
- Returns:
int
- Return type:
estimated embedding dimension
- att.embedding.validate_embedding(cloud, expected_dim=None, condition_threshold=10000.0)[source]¶
Check embedding quality via SVD of the centered point cloud matrix.
The condition number is σ_max / σ_min. High values mean some embedding dimensions are near-linear combinations of others — the manifold is collapsed along those directions.
- Parameters:
cloud ((n_points, dimension) array)
expected_dim (expected intrinsic dimension (informational only))
condition_threshold (condition number above which embedding is flagged) – as degenerate. Default 1e4 calibrated on coupled Rössler-Lorenz.
- Returns:
condition_number: float singular_values: np.ndarray effective_rank: int (singular values > 1e-3 * σ_max) degenerate: bool warning: str | None
- Return type:
dict with keys
- att.embedding.svd_embedding(X, delay, dimension, n_components=None)[source]¶
SVD-projected delay embedding for noise reduction.
Constructs the delay matrix then projects onto the top n_components principal components.
- Parameters:
X ((n_samples,) 1D time series)
delay (time delay)
dimension (embedding dimension)
n_components (number of SVD components to keep (default: dimension))
- Return type:
(n_valid, n_components) projected point cloud
- exception att.embedding.EmbeddingDegeneracyWarning[source]¶
Bases:
UserWarningRaised when an embedding has a dangerously high condition number.
- class att.embedding.TakensEmbedder(delay='auto', dimension='auto')[source]¶
Bases:
objectReconstruct a phase-space attractor from a scalar time series.
- Parameters:
- fit(X)[source]¶
Estimate parameters from data. Stores .delay_ and .dimension_.
- Parameters:
X (ndarray)
- Return type:
- class att.embedding.JointEmbedder(delays='auto', dimensions='auto')[source]¶
Bases:
objectConstruct joint delay embeddings with per-channel delay estimation.
Using “auto” is strongly recommended for systems with different timescales.
- Parameters:
- fit(channels)[source]¶
Estimate per-channel parameters. Stores .delays_ and .dimensions_.
- Parameters:
- Return type:
- transform(channels)[source]¶
Construct joint delay vectors by concatenating per-channel embeddings.
Input: list of 1D arrays, each (n_samples,) Output: (n_valid_samples, sum(dimensions))
Time-delay estimation via Average Mutual Information.
- att.embedding.delay.estimate_delay(X, method='ami', max_lag=100)[source]¶
Estimate optimal time delay for Takens embedding.
Uses first minimum of Average Mutual Information (Fraser & Swinney, 1986).
- Parameters:
X ((n_samples,) array)
method ("ami" (only option currently))
max_lag (maximum lag to consider)
- Returns:
int
- Return type:
optimal delay (first AMI minimum, or max_lag if no minimum found)
Embedding dimension estimation via False Nearest Neighbors.
- att.embedding.dimension.estimate_dimension(X, delay, method='fnn', max_dim=10, threshold=0.01, rtol=15.0, atol=2.0)[source]¶
Estimate minimal embedding dimension via False Nearest Neighbors.
- Parameters:
X ((n_samples,) 1D time series)
delay (time delay for embedding)
method ("fnn" (only option currently))
max_dim (maximum dimension to test)
threshold (FNN fraction below which to stop (default 0.01 = 1%))
rtol (relative distance threshold for FNN criterion 1)
atol (absolute distance threshold for FNN criterion 2)
- Returns:
int
- Return type:
estimated embedding dimension
Embedding quality validation via condition number analysis.
- class att.embedding.validation.EmbeddingDegeneracyWarning[source]¶
Bases:
UserWarningRaised when an embedding has a dangerously high condition number.
- exception att.embedding.validation.EmbeddingDegeneracyWarning[source]¶
Bases:
UserWarningRaised when an embedding has a dangerously high condition number.
- att.embedding.validation.validate_embedding(cloud, expected_dim=None, condition_threshold=10000.0)[source]¶
Check embedding quality via SVD of the centered point cloud matrix.
The condition number is σ_max / σ_min. High values mean some embedding dimensions are near-linear combinations of others — the manifold is collapsed along those directions.
- Parameters:
cloud ((n_points, dimension) array)
expected_dim (expected intrinsic dimension (informational only))
condition_threshold (condition number above which embedding is flagged) – as degenerate. Default 1e4 calibrated on coupled Rössler-Lorenz.
- Returns:
condition_number: float singular_values: np.ndarray effective_rank: int (singular values > 1e-3 * σ_max) degenerate: bool warning: str | None
- Return type:
dict with keys
- att.embedding.validation.svd_embedding(X, delay, dimension, n_components=None)[source]¶
SVD-projected delay embedding for noise reduction.
Constructs the delay matrix then projects onto the top n_components principal components.
- Parameters:
X ((n_samples,) 1D time series)
delay (time delay)
dimension (embedding dimension)
n_components (number of SVD components to keep (default: dimension))
- Return type:
(n_valid, n_components) projected point cloud