epyr.signalprocessing

Signal Processing Module for EPR Time-Domain Data

This module provides frequency-analysis tools for time-dependent EPR signals, with optimisations for pulse-EPR experiments: Rabi oscillations, DEER, HYSCORE, and related techniques.

Key Features

  • 1D FFT analysis with automatic DC offset removal and windowing

  • 2D FFT analysis with two processing modes: - Row-by-row 1D FFT for 2D Rabi oscillations and similar experiments - Full 2D FFT for HYSCORE-type measurements with cross-peak analysis

  • Advanced apodization windows for spectral leakage reduction

  • Power spectral density analysis using Welch and periodogram methods

  • Time-frequency analysis with spectrograms

  • Automatic time unit detection (ns, μs, ms, s)

  • Zero padding for enhanced frequency resolution

Main Functions

analyze_frequencies : 1D FFT-based frequency analysis with four-panel plot analyze_frequencies_2d : 2D FFT analysis with row-by-row or full 2D processing modes power_spectrum : Power spectral density using Welch or periodogram methods spectrogram_analysis : Time-frequency analysis for evolving spectral content apowin : Apodization window generation with multiple window types

Examples

Basic 1D frequency analysis of Rabi data:

from epyr import eprload
from epyr.signalprocessing import analyze_frequencies

# Load time-domain EPR data
time, signal, params, _ = eprload('rabi_data.DTA')

# Analyze frequencies with DC removal and Hann window
result = analyze_frequencies(time, signal, window='hann',
                           remove_dc=True, plot=True)

print(f"Dominant frequency: {result['dominant_frequencies'][0]:.3f} MHz")

2D frequency analysis - Row-by-row mode (2D Rabi oscillations):

from epyr.signalprocessing import analyze_frequencies_2d

# Load 2D Rabi data
x_2d, y_2d, params, _ = eprload('rabi_2d.DTA')

# Process each trace independently with 1D FFT
result = analyze_frequencies_2d(x_2d[0], y_2d, mode='row_by_row',
                               window='hann', plot=True)

print(f"Power spectrum shape: {result['power_spectrum'].shape}")

2D frequency analysis - Full 2D FFT (HYSCORE):

from epyr.signalprocessing import analyze_frequencies_2d

# Load HYSCORE data
x_hyscore, y_hyscore, params, _ = eprload('hyscore.DTA')

# Apply full 2D FFT for cross-peak analysis
result = analyze_frequencies_2d((x_hyscore[0], x_hyscore[1]), y_hyscore,
                               mode='full_2d', window='hann', plot=True)

print(f"2D power spectrum shape: {result['power_spectrum'].shape}")

Advanced power spectrum analysis:

from epyr.signalprocessing import power_spectrum

# Welch method for noise reduction
psd_result = power_spectrum(time, signal, method='welch',
                          window='hann', plot=True)

Notes

This module is specifically designed for EPR time-domain signal analysis and includes optimizations for typical EPR data characteristics including proper handling of complex signals, automatic unit detection, and spectroscopic conventions.

The 2D FFT capabilities enable analysis of: - 2D Rabi oscillations (row-by-row mode) - HYSCORE spectra (full 2D mode) - Other multi-dimensional time-domain EPR experiments

epyr.signalprocessing.apowin(window_type, n_points, alpha=None, half_window=None)[source]

Generate apodization windows for signal processing.

Apodization windows are used to reduce spectral leakage and improve signal-to-noise ratio in Fourier transform spectroscopy.

Parameters:

window_typestr

Window type: - ‘hamming’ or ‘ham’: Hamming window - ‘hann’ or ‘han’: Hann (Hanning) window - ‘blackman’ or ‘bla’: Blackman window - ‘bartlett’ or ‘bar’: Bartlett (triangular) window - ‘connes’ or ‘con’: Connes window - ‘cosine’ or ‘cos’: Cosine window - ‘welch’ or ‘wel’: Welch window - ‘kaiser’ or ‘kai’: Kaiser window (needs alpha) - ‘gaussian’ or ‘gau’: Gaussian window (needs alpha) - ‘exponential’ or ‘exp’: Exponential window (needs alpha)

n_pointsint

Number of points in the window

alphafloat, optional

Shape parameter for Kaiser, Gaussian, and Exponential windows

half_windowstr, optional

Generate half window: ‘left’ (-1 to 0) or ‘right’ (0 to 1)

Returns:

array

Normalized window values (peak = 1)

Examples:

>>> # Hamming window
>>> w = apowin('hamming', 256)
>>> # Kaiser window with beta=6
>>> w_kaiser = apowin('kaiser', 256, alpha=6)
>>> # Half Hann window (right side)
>>> w_half = apowin('hann', 128, half_window='right')
epyr.signalprocessing.window_comparison(n_points=256)[source]

Compare different window types side by side.

Parameters:

n_pointsint

Number of points for each window

epyr.signalprocessing.frequency_response_demo(n_points=256)[source]

Show frequency response characteristics of different windows.

epyr.signalprocessing.apply_window_demo()[source]

Demonstrate windowing effect on a test signal

epyr.signalprocessing.analyze_frequencies(time_data, signal_data, window='hann', window_alpha=None, zero_padding=2, remove_dc=True, plot=True, freq_range=None, **plot_kwargs)[source]

FFT-based frequency analysis of time-domain EPR signals.

This function performs clean FFT analysis to identify frequency components in time-dependent EPR signals, with proper DC offset removal.

Parameters:

time_datanp.ndarray

Time axis data (in ns, μs, or s)

signal_datanp.ndarray

EPR signal intensity vs time

windowstr or None, optional

Apodization window type (‘hann’, ‘hamming’, ‘blackman’, ‘kaiser’, None) Default: ‘hann’

window_alphafloat, optional

Alpha parameter for Kaiser, Gaussian windows (default: 6 for Kaiser)

zero_paddingint, optional

Zero padding factor (2 = double length, 4 = quadruple, etc.) Default: 2

remove_dcbool, optional

Remove DC offset before analysis (recommended: True)

plotbool, optional

Generate analysis plots. Default: True

freq_rangetuple of float, optional

Frequency range (min, max) to display in plots

Returns:

dict

Analysis results containing: - ‘frequencies’: Frequency axis in appropriate units - ‘power_spectrum’: Power spectral density (normalized) - ‘phase_spectrum’: Phase spectrum - ‘dominant_frequencies’: List of peak frequencies - ‘time_data’: Original time data - ‘processed_signal’: Signal after DC removal and windowing - ‘sampling_rate’: Sampling rate in Hz - ‘time_unit’: Detected time unit - ‘freq_unit’: Frequency unit

Examples:

>>> from epyr import eprload
>>> from epyr.signalprocessing import analyze_frequencies
>>>
>>> # Load Rabi data
>>> time, signal, params, _ = eprload('rabi_data.DTA')
>>> result = analyze_frequencies(time, signal, window='hann', plot=True)
>>> print(f"Dominant frequency: {result['dominant_frequencies'][0]:.3f} MHz")
Parameters:
Return type:

Dict

epyr.signalprocessing.analyze_frequencies_2d(time_data, signal_data, mode='row_by_row', window='hann', window_alpha=None, zero_padding=2, remove_dc=True, axis=1, plot_result=False, freq_range=None, **plot_kwargs)[source]

FFT-based frequency analysis of 2D time-domain EPR signals.

This function handles 2D EPR data with two processing modes: 1. Row-by-row 1D FFT: Process each row/column independently (e.g., 2D Rabi) 2. Full 2D FFT: Process both dimensions together (e.g., HYSCORE)

Parameters:

time_datanp.ndarray or tuple of np.ndarray

Time axis data. Can be: - Single 1D array: time axis for the FFT dimension - Tuple of two 1D arrays: (time_axis1, time_axis2) for 2D FFT

signal_datanp.ndarray

2D EPR signal intensity array (shape: n_traces x n_points)

modestr, optional

Processing mode: - ‘row_by_row’: Apply 1D FFT to each row/column independently - ‘full_2d’: Apply 2D FFT to entire dataset (HYSCORE-type) Default: ‘row_by_row’

windowstr or None, optional

Apodization window type. Default: ‘hann’

window_alphafloat, optional

Alpha parameter for Kaiser, Gaussian windows

zero_paddingint, optional

Zero padding factor. Default: 2

remove_dcbool, optional

Remove DC offset before analysis. Default: True

axisint, optional

Axis to process for row_by_row mode (0=columns, 1=rows). Default: 1

plot_resultbool, optional

Generate analysis plots. Default: False

freq_rangetuple of float, optional

Frequency range (min, max) to display in plots

Returns:

For mode=’row_by_row’:
fqnp.ndarray

Frequency axis (1D array)

axis2np.ndarray

Secondary axis (field, angle, trace index, etc.)

spectrumnp.ndarray

2D FFT spectrum magnitude (n_traces x n_frequencies)

infodict

Analysis information (units, sampling_rate, mode, etc.)

For mode=’full_2d’:
fq1np.ndarray

Frequency axis 1 (1D array)

fq2np.ndarray

Frequency axis 2 (1D array)

spectrumnp.ndarray

2D FFT spectrum magnitude (n_freq1 x n_freq2)

infodict

Analysis information

Examples:

>>> # Row-by-row 1D FFT (2D Rabi oscillations)
>>> x_2d, y_2d, params, _ = eprload('rabi_2d.DTA')
>>> fq, axis2, spectrum, info = analyze_frequencies_2d(
...     x_2d[0], y_2d, mode='row_by_row', plot_result=False)
>>> # Full 2D FFT (HYSCORE)
>>> x_hyscore, y_hyscore, params, _ = eprload('hyscore.DTA')
>>> fq1, fq2, spectrum_2d, info = analyze_frequencies_2d(
...     (x_hyscore[0], x_hyscore[1]), y_hyscore,
...     mode='full_2d', plot_result=True)
Parameters:
Return type:

Tuple

epyr.signalprocessing.power_spectrum(time_data, signal_data, method='welch', window='hann', nperseg=None, overlap=0.5, remove_dc=True, plot=True)[source]

Calculate power spectral density using Welch or periodogram methods.

Parameters:

time_datanp.ndarray

Time axis data

signal_datanp.ndarray

Signal data

methodstr

Method: ‘welch’ or ‘periodogram’

windowstr

Window function for Welch method

npersegint, optional

Length of each segment for Welch method

overlapfloat

Overlap fraction for Welch method (0-1)

remove_dcbool

Remove DC offset before analysis

plotbool

Generate plots

Returns:

dict

Results with frequencies and power spectrum

Parameters:
Return type:

Dict

epyr.signalprocessing.spectrogram_analysis(time_data, signal_data, window='hann', nperseg=None, overlap=0.8, remove_dc=True, plot=True)[source]

Time-frequency analysis using spectrogram.

Parameters:

time_datanp.ndarray

Time axis data

signal_datanp.ndarray

Signal data

windowstr

Window function

npersegint, optional

Length of each segment

overlapfloat

Overlap fraction (0-1)

remove_dcbool

Remove DC offset

plotbool

Generate spectrogram plot

Returns:

dict

Results with time axis, frequencies, and spectrogram

Parameters:
Return type:

Dict

Modules

apowin(window_type, n_points[, alpha, ...])

Generate apodization windows for signal processing.

frequency_analysis

Frequency Analysis for Time-Domain EPR Signals