impedance_agent.fitters package
Submodules
impedance_agent.fitters.drt module
- class DRTFitter(zexp_re, zexp_im, omg, lam_t0, lam_pg0, lower_bounds, upper_bounds, mode='real')[source]
Bases:
object
Distribution of Relaxation Times (DRT) analysis for electrochemical impedance data.
This class implements the combined Tikhonov regularization and projected gradient method for calculating DRT from impedance data as described in: “PEM fuel cell distribution of relaxation times: a method for the calculation and behavior of an oxygen transport peak”
The fundamental DRT equation relating impedance to the distribution function is:
\[Z(\omega) - R_\infty = R_{pol}\int_0^\infty \frac{\gamma(\tau)d\tau}{1 + i\omega\tau}\]With normalization condition:
\[\int_0^\infty \gamma(\tau)d\tau = 1\]- Parameters:
zexp_re (np.ndarray) – Real part of experimental impedance
zexp_im (np.ndarray) – Imaginary part of experimental impedance
omg (np.ndarray) – Angular frequencies
lam_t0 (float) – Initial Tikhonov regularization parameter
lam_pg0 (float) – Initial projected gradient parameter
lower_bounds (np.ndarray) – Lower bounds for regularization parameters
upper_bounds (np.ndarray) – Upper bounds for regularization parameters
mode (str, optional) – ‘real’ or ‘imag’ for using real or imaginary part of impedance (default: ‘real’)
- compute_normalized_residuals(Z_fit: ndarray) Tuple[ndarray, ndarray] [source]
Compute normalized residuals using impedance modulus for normalization.
\[ \begin{align}\begin{aligned}r_{real} = \frac{Z_{fit,real} - Z_{exp,real}}{\sqrt{Z_{exp,real}^2 + Z_{exp,imag}^2}}\\r_{imag} = \frac{Z_{fit,imag} - Z_{exp,imag}}{\sqrt{Z_{exp,real}^2 + Z_{exp,imag}^2}}\end{aligned}\end{align} \]- Parameters:
Z_fit (np.ndarray) – Complex array containing the fitted impedance values
- Returns:
Tuple containing: - residuals_real: Normalized residuals of the real component - residuals_imag: Normalized residuals of the imaginary component
- Return type:
Tuple[np.ndarray, np.ndarray]
- tikh_solver(lam_t, a_mat_t_a, b_rhs, id_matrix)[source]
Solve the Tikhonov-regularized equation using least squares.
\[(A^TA + \lambda_T I)\vec{\gamma} = A^T\vec{b}\]- Parameters:
lam_t (float) – Tikhonov regularization parameter
a_mat_t_a (np.ndarray) – A^T A matrix
b_rhs (np.ndarray) – Right-hand side vector
id_matrix (np.ndarray) – Identity matrix
- Returns:
Solution vector γ
- Return type:
np.ndarray
- objective_function(g_vector, lhs_matrix, b_rhs)[source]
Objective function for projected gradient optimization.
\[f(\vec{\gamma}) = \|(A^TA + \lambda_T I)\vec{\gamma} - A^T\vec{b}\|^2\]- Parameters:
g_vector (np.ndarray) – Current γ vector
lhs_matrix (np.ndarray) – Left-hand side matrix
b_rhs (np.ndarray) – Right-hand side vector
- Returns:
Objective function value
- Return type:
- find_lambda()[source]
Scans over lambda_T range to find optimal regularization parameters.
- Returns:
resid: Residual norms
solnorm: Solution norms
lam_t_arr: Tikhonov parameter values
lam_pg_arr: Projected gradient parameter values
- Return type:
Tuple[np.ndarray, np.ndarray, np.ndarray, np.ndarray]
- residual_norm(g_vector)[source]
Calculate the norm of the residual using (A^T A)g - A^T b.
- Parameters:
g_vector (np.ndarray) – Current γ vector
- Returns:
Residual norm
- Return type:
- encode(p, lb, ub)[source]
Converts external parameters to internal parameters in log10 scale.
\[p_{int} = \log_{10}\left(\frac{p - lb}{1 - p/ub}\right)\]
- decode(p, lb, ub)[source]
Converts internal parameters to external parameters.
\[p_{ext} = \frac{lb + 10^p}{1 + 10^p/ub}\]
- jacobian_lsq(pvec, lhs_matrix, a_mat_t_a, b_rhs, d_ln_tau, id_matrix, lb, ub)[source]
Compute the Jacobian of the Tikhonov residual function.
- pg_solver(lamvec, lhs_matrix, a_mat_t_a, b_rhs, d_ln_tau, id_matrix)[source]
Projected gradient solver using Tikhonov solution as initial guess.
Implements the iterations:
\[\vec{\gamma}^{k+1} = P_+[\vec{\gamma}^k - \lambda_{pg}((A^TA + \lambda_TI)\vec{\gamma}^k - A^T\vec{b})]\]where P_+ is the projection onto non-negative orthant.
- tikh_residual(lamvec_log, lhs_matrix, a_mat_t_a, b_rhs, d_ln_tau, id_matrix, lb, ub)[source]
Compute Tikhonov residual for given regularization parameters.
- rpol_peaks(g_vector)[source]
Find peaks in the DRT spectrum and calculate their parameters.
- Parameters:
g_vector (np.ndarray) – DRT solution vector
- Returns:
Array containing peak parameters (frequencies and polarizations)
- Return type:
np.ndarray
- z_model_imre(g_vector)[source]
Calculate the model impedance from the DRT solution.
\[Z(\omega) = R_{pol}\int_0^\infty \frac{\gamma(\tau)}{1 + i\omega\tau}d\tau\]
- fit() DRTResult | None [source]
Perform the complete DRT fitting process.
- Returns:
Complete DRT analysis results or None if fitting fails
- Return type:
Optional[DRTResult]
Notes
Finds optimal regularization parameters
Performs Tikhonov regularization with projected gradient
Analyzes peaks in the DRT spectrum
Calculates fitted impedance and residuals
impedance_agent.fitters.ecm module
- class ECMFitter(model_func, p0, freq, impedance_data: ImpedanceData, lb, ub, param_info, weighting='modulus')[source]
Bases:
object
Equivalent Circuit Model (ECM) fitting for electrochemical impedance data.
This class implements least squares optimization with bounded parameters for fitting equivalent circuit models to impedance data. The fitting process uses weighted residuals and supports different weighting schemes.
The fundamental equation for the weighted sum of squared residuals is:
\[WRSS = \sum_{i=1}^N \frac{(Z_{\mathrm{exp},i} - Z_{\mathrm{model},i})^2}{\sigma_i^2}\]Supported weighting schemes:
\[\begin{split}\sigma_i = \begin{cases} 1 & \text{for unit weighting} \\ |Z_{\mathrm{exp},i}| & \text{for proportional weighting} \\ \sqrt{(\mathrm{Re}(Z_{\mathrm{exp},i}))^2 + (\mathrm{Im}(Z_{\mathrm{exp},i}))^2} & \text{for modulus weighting} \end{cases}\end{split}\]- Parameters:
model_func (callable) – Function that takes parameters and frequencies and returns impedance
p0 (array_like) – Initial parameter values
freq (array_like) – Frequency values
impedance_data (ImpedanceData) – Object containing experimental impedance data
lb (array_like) – Lower bounds for parameters
ub (array_like) – Upper bounds for parameters
param_info (list) – List of dictionaries containing parameter information
weighting (str or array_like, optional) – Weighting scheme (‘unit’, ‘proportional’, ‘modulus’) or custom weights (default: ‘modulus’)
- __init__(model_func, p0, freq, impedance_data: ImpedanceData, lb, ub, param_info, weighting='modulus')[source]
Initialize ECM fitter with model and data.
- Parameters:
model_func (callable) – Function that takes parameters and frequencies and returns impedance
p0 (array_like) – Initial parameter values
freq (array_like) – Frequency values
impedance_data (ImpedanceData) – Object containing experimental impedance data
lb (array_like) – Lower bounds for parameters
ub (array_like) – Upper bounds for parameters
param_info (list) – List of dictionaries containing parameter information
weighting (str or array_like, optional) – Weighting scheme (‘unit’, ‘proportional’, ‘modulus’) or custom weights (default: ‘modulus’)
- compute_normalized_residuals(Z_fit: ndarray) Tuple[ndarray, ndarray] [source]
Compute normalized residuals using impedance modulus for normalization.
\[ \begin{align}\begin{aligned}r_{real} = \frac{Z_{fit,real} - Z_{exp,real}}{\sqrt{Z_{exp,real}^2 + Z_{exp,imag}^2}}\\r_{imag} = \frac{Z_{fit,imag} - Z_{exp,imag}}{\sqrt{Z_{exp,real}^2 + Z_{exp,imag}^2}}\end{aligned}\end{align} \]- Parameters:
Z_fit (np.ndarray) – Complex array containing the fitted impedance values
- Returns:
Tuple containing: - residuals_real: Normalized residuals of the real component - residuals_imag: Normalized residuals of the imaginary component
- Return type:
Tuple[np.ndarray, np.ndarray]
- encode(p: Array) Array [source]
Convert external parameters to internal parameters using log-transform.
\[p_{int} = \log_{10}\left(\frac{p - lb}{1 - p/ub}\right)\]- Parameters:
p (jnp.ndarray) – External parameters
- Returns:
Internal parameters
- Return type:
jnp.ndarray
- decode(p: Array) Array [source]
Convert internal parameters to external parameters.
\[p_{ext} = \frac{lb + 10^p}{1 + 10^p/ub}\]- Parameters:
p (jnp.ndarray) – Internal parameters
- Returns:
External parameters
- Return type:
jnp.ndarray
- obj_fun(p_log: Array) float [source]
Calculate weighted objective function.
\[WRSS = \sum_{i=1}^N \frac{(Z_{exp,i} - Z_{model,i})^2}{\sigma_i^2}\]- Parameters:
p_log (jnp.ndarray) – Parameters in log space
- Returns:
Weighted residual sum of squares
- Return type:
- compute_aic(wrss: float) float [source]
Compute Akaike Information Criterion (AIC) for model selection.
For unit weighting:
\[\mathrm{AIC} = 2N\ln(2\pi) - 2N\ln(2N) + 2N + 2N\ln(\mathrm{WRSS}) + 2k\]For modulus/proportional weighting:
\[\mathrm{AIC} = 2N\ln(2\pi) - 2N\ln(2N) + 2N - \sum\ln(w_i) + 2N\ln(\mathrm{WRSS}) + 2(k+1)\]For sigma weighting:
\[\mathrm{AIC} = 2N\ln(2\pi) + \sum\ln(\sigma_i^2) + \mathrm{WRSS} + 2k\]
- calculate_fitted_impedance(parameters: Array) Array [source]
Calculate fitted impedance values from parameters.
- Parameters:
parameters (jnp.ndarray) – Model parameters
- Returns:
Complex array of fitted impedance values
- Return type:
jnp.ndarray
- compute_fit_quality_metrics(Z_fit: ndarray) FitQualityMetrics [source]
Compute fit quality using vector difference and path deviation metrics.
The vector difference analysis quantifies the point-by-point agreement:
\[\mathrm{VD} = \frac{1}{N}\sum_{i=1}^N \frac{|Z_{\mathrm{fit},i} - Z_{\mathrm{exp},i}|}{|Z_{\mathrm{exp},i}|}\]The path deviation analysis quantifies trajectory agreement:
\[\mathrm{PD} = \frac{1}{N-1}\sum_{i=1}^{N-1} \left|\frac{\Delta Z_{\mathrm{fit},i}}{|\Delta Z_{\mathrm{fit},i}|} - \frac{\Delta Z_{\mathrm{exp},i}}{|\Delta Z_{\mathrm{exp},i}|}\right|\]- Parameters:
Z_fit (np.ndarray) – Complex fitted impedance array
- Returns:
Computed quality metrics including vector difference, path deviation, and overall quality assessment
- Return type:
- fit() FitResult | None [source]
Perform impedance fitting on the data.
The fitting process involves: 1. Parameter optimization using BFGS algorithm 2. Uncertainty calculation via QR decomposition 3. Computation of fit quality metrics 4. Calculation of correlation matrix 5. Model selection metrics (AIC)
- Returns:
Complete fitting results including: - Optimized parameters and uncertainties - Correlation matrix - Goodness-of-fit metrics (χ², AIC, WRMS) - Fitted impedance values - Fit quality assessment Returns None if fitting fails
- Return type:
Optional[FitResult]
impedance_agent.fitters.linkk module
- class LinKKFitter(data: ImpedanceData)[source]
Bases:
object
- __init__(data: ImpedanceData)[source]