Homodyne Scattering Analysis Package
A high-performance Python package for analyzing homodyne scattering in X-ray Photon Correlation Spectroscopy (XPCS) under nonequilibrium conditions. Implements the theoretical framework from He et al. PNAS 2024 for characterizing transport properties in flowing soft matter systems.
Overview
This package analyzes time-dependent intensity correlation functions c₂(φ,t₁,t₂) in complex fluids under nonequilibrium conditions, capturing the interplay between Brownian diffusion and advective shear flow. The implementation provides:
Three analysis modes: Static Isotropic (3 params), Static Anisotropic (3 params), Laminar Flow (7 params)
Multiple optimization methods: Classical (Nelder-Mead, Iterative Gurobi with Trust Regions), Robust (Wasserstein DRO, Scenario-based, Ellipsoidal)
High performance: Numba JIT compilation with 3-5x speedup, vectorized NumPy operations, comprehensive performance monitoring
Scientific accuracy: Automatic c₂ = offset + contrast × c₁ fitting for proper chi-squared calculations
Quick Start
Installation:
pip install homodyne-analysis[all]
Python API:
import numpy as np
import json
from homodyne.analysis.core import HomodyneAnalysisCore
from homodyne.optimization.classical import ClassicalOptimizer
from homodyne.data.xpcs_loader import load_xpcs_data
# Load configuration
with open("config.json", 'r') as f:
config = json.load(f)
# Initialize analysis core
core = HomodyneAnalysisCore(config)
# Load experimental data
phi_angles = np.array([0, 36, 72, 108, 144])
c2_data = load_xpcs_data(
data_path=config['experimental_data']['data_folder_path'],
phi_angles=phi_angles,
n_angles=len(phi_angles)
)
# Run optimization
optimizer = ClassicalOptimizer(core, config)
params, results = optimizer.run_classical_optimization_optimized(
phi_angles=phi_angles,
c2_experimental=c2_data
)
print(f"D₀ = {params[0]:.3e} Ų/s")
print(f"χ² = {results.chi_squared:.6e}")
Command Line Interface:
# Configuration generator
homodyne-config --mode static_isotropic --sample protein_01
homodyne-config --mode laminar_flow --sample microgel
# Main analysis command
homodyne # Default classical method
homodyne --method robust # Robust optimization only
homodyne --method all --verbose # All methods with debug logging
# Analysis mode control
homodyne --static-isotropic # Force 3-parameter isotropic mode
homodyne --static-anisotropic # Force 3-parameter anisotropic mode
homodyne --laminar-flow # Force 7-parameter flow mode
# Data visualization
homodyne --plot-experimental-data # Validate experimental data
homodyne --plot-simulated-data # Plot theoretical correlations
homodyne --plot-simulated-data --contrast 1.5 --offset 0.1 --phi-angles "0,45,90,135"
# Configuration and output
homodyne --config my_config.json --output-dir ./results --verbose
homodyne --quiet # File logging only, no console output
What’s New in v1.0.0
Critical Bug Fixes:
Frame Counting Convention - Fixed 1-based inclusive counting to 0-based Python slicing conversion
Formula:
time_length = end_frame - start_frame + 1(inclusive)Resolves NaN chi-squared values and cache dimension mismatches
Applied across all 9 modules: analysis core, data loader, optimizers, CLI
Conditional Angle Subsampling - Preserve angular information when
n_angles < 4Prevents loss of critical angular information (e.g., 2 angles → 1 angle)
Time subsampling still applied for performance (~16x reduction)
Implemented in both classical and robust optimizers
Memory Optimization - Increased ellipsoidal optimization memory limit to 90%
Handles large datasets with 8M+ data points without overflow
Fixed stacked decorator issue in robust optimization
Performance Improvements:
Numba JIT Compilation - 3-5x speedup for core calculations
Vectorized Operations - Optimized NumPy array processing throughout
Smart Caching - Intelligent data caching with automatic dimension validation
Analysis Modes
Mode |
Parameters |
Use Case |
Speed |
Command |
|---|---|---|---|---|
Static Isotropic |
3 |
Fastest, isotropic systems |
⭐⭐⭐ |
|
Static Anisotropic |
3 |
Static with angular dependencies |
⭐⭐ |
|
Laminar Flow |
7 |
Flow & shear analysis |
⭐ |
|
Key Features
- Multiple Analysis Modes
Static Isotropic (3 parameters), Static Anisotropic (3 parameters), and Laminar Flow (7 parameters)
- High Performance
Numba JIT compilation, smart angle filtering, and optimized computational kernels
- Scientific Accuracy
Automatic c₂ = offset + contrast × c₁ fitting for accurate chi-squared calculations
Multiple Optimization Methods
- Security and Code Quality
Comprehensive security scanning with Bandit, dependency vulnerability checking with pip-audit, and automated code quality tools
- Comprehensive Validation
Experimental data validation plots and quality control
Visualization Tools
- Performance Monitoring
Comprehensive performance testing, regression detection, and automated benchmarking
User Guide
- Installation Guide
- Quick Start Guide
- Analysis Modes
- Configuration Guide
- Data Visualization
- ML Training Data - Fixed and Ready! 🎉
- Examples
- Example 1: Basic Isotropic Analysis
- Example 2: Flow Analysis with Robust Methods
- Example 3: Performance-Optimized Analysis
- Example 4: Batch Processing Multiple Samples
- Example 5: Progressive Analysis Workflow
- Common Patterns
- Output Directory Structure
- Diagnostic Summary Visualizations
- Common Output Structure for All 5 Classical Methods
- Available Optimization Methods
- Example 6: Logging Control for Different Scenarios
- Example 7: Performance Monitoring and Optimization
- Next Steps
API Reference
Developer Guide
- Contributing Guide
- Testing Guide
- Performance Optimization
- Performance Guide
- Performance Overview (v0.6.5+)
- Method Performance Comparison
- Performance Optimization Strategies
- Performance Monitoring
- Performance Testing
- Environment Optimization
- Troubleshooting Performance Issues
- Best Practices
- Code Quality and Maintenance
- Performance Overview
- Optimization Strategies
- Memory Optimization
- CPU Optimization
- Algorithm Optimization
- Performance Benchmarks
- Profiling Tools
- Performance Best Practices
- Troubleshooting Performance Issues
- Package Architecture
- Troubleshooting Guide
Theoretical Background
The package implements three key equations describing correlation functions in nonequilibrium laminar flow systems:
- Equation 13 - Full Nonequilibrium Laminar Flow:
c₂(q⃗, t₁, t₂) = 1 + β[e^(-q²∫J(t)dt)] × sinc²[1/(2π) qh ∫γ̇(t)cos(φ(t))dt]
- Equation S-75 - Equilibrium Under Constant Shear:
c₂(q⃗, t₁, t₂) = 1 + β[e^(-6q²D(t₂-t₁))] sinc²[1/(2π) qh cos(φ)γ̇(t₂-t₁)]
- Equation S-76 - One-time Correlation (Siegert Relation):
c₂(q⃗, τ) = 1 + β[e^(-6q²Dτ)] sinc²[1/(2π) qh cos(φ)γ̇τ]
Key Parameters:
q⃗: scattering wavevector [Å⁻¹]
h: gap between stator and rotor [Å]
φ(t): angle between shear/flow direction and q⃗ [degrees]
γ̇(t): time-dependent shear rate [s⁻¹]
D(t): time-dependent diffusion coefficient [Ų/s]
β: contrast parameter [dimensionless]
Citation
If you use this package in your research, please cite:
@article{he2024transport,
title={Transport coefficient approach for characterizing nonequilibrium dynamics in soft matter},
author={He, Hongrui and Liang, Hao and Chu, Miaoqi and Jiang, Zhang and de Pablo, Juan J and Tirrell, Matthew V and Narayanan, Suresh and Chen, Wei},
journal={Proceedings of the National Academy of Sciences},
volume={121},
number={31},
pages={e2401162121},
year={2024},
publisher={National Academy of Sciences},
doi={10.1073/pnas.2401162121}
}
Support
Documentation: https://homodyne.readthedocs.io/
Source Code: https://github.com/imewei/homodyne
License: MIT License