impedance.py¶
impedance.py
is a Python module for working with impedance data.
This project started at the 2018 Electrochemical Society (ECS) Hack Week and has grown from there.
Using a scikit-learn-like API, we hope to make visualizing, fitting, and analyzing impedance spectra more intuitive and reproducible.
Note
impedance.py
is currently in a beta phase and new features are rapidly being added.
If you have a feature request or find a bug, please feel free to file an issue or, better yet, make the code improvements and submit a pull request! The goal is to build an open-source tool that the entire impedance community can improve and use!
User Installation¶
The easiest way to install impedance.py
is from PyPI using pip:
pip install impedance
Dependencies¶
impedance.py requires:
- Python (>=3.5)
- SciPy (>=1.0)
- NumPy (>=1.14)
- Matplotlib (>=3.0)
Several example notebooks are provided in the examples/
directory.
Opening these will also require Jupyter Notebook or Jupyter Lab.
Examples and Documentation¶
Getting started with impedance.py contains a detailed walk
through of how to get started from scratch. If you’re already familiar with
Jupyter/Python, several examples can be found in the examples/
directory
(the fitting_tutorial.ipynb
is a great place to start) and the documentation can be found at
impedancepy.readthedocs.io.
Getting started with impedance.py
¶
impedance.py
is a Python package for analyzing electrochemical impedance
spectroscopy (EIS) data. The following steps will show you how to get started
analyzing your own data using impedance.py
in a Jupyter notebook.
Hint
If you get stuck or believe you have found a bug, please feel free to open an issue on GitHub.
Step 1: Installation¶
If you already are familiar with managing Python packages, feel free to skip straight to Installing packages. Otherwise, what follows is a quick introduction to the Python package ecosystem:
Installing Miniconda¶
One of the easiest ways to get started with Python is using Miniconda. Installation instructions for your OS can be found at https://conda.io/miniconda.html.
After you have installed conda, you can run the following commands in your Terminal/command prompt/Git BASH to update and test your installation:
- Update conda’s listing of packages for your system:
conda update conda
- Test your installation:
conda list
For a successful installation, a list of installed packages appears.
- Test that Python 3 is your default Python:
python -V
You should see something like Python 3.x.x :: Anaconda, Inc.
You can interact with Python at this point, by simply typing python
.
Setting up a conda environment¶
(Optional) It is recommended that you use virtual environments to keep track of the packages you’ve installed for a particular project. Much more info on how conda makes this straightforward is given here.
We will start by creating an environment called impedance-analysis
which contains all the Python base distribution:
conda create -n impedance-analysis python=3
After conda creates this environment, we need to activate it before we can install anything into it by using:
conda activate impedance-analysis
We’ve now activated our conda environment and are ready to install
impedance.py
!
Installing packages¶
The easiest way to install impedance.py
and it’s dependencies
(scipy
, numpy
, and matplotlib
) is from
PyPI using pip:
pip install impedance
For this example we will also need Jupyter Lab which we can install with:
conda install jupyter jupyterlab
We’ve now got everything in place to start analyzing our EIS data!
Note
The next time you want to use this same environment, all you have to do is
open your terminal and type conda activate impedance-analysis
.
Open Jupyter Lab¶
(Optional) Create a directory in your documents folder for this example:
mkdir ~/Documents/impedance-example
cd ~/Documents/impedance-example
Next, we will launch an instance of Jupyter Lab:
jupyter lab
which should open a new tab in your browser. A tutorial on Jupyter Lab from the Electrochemical Society HackWeek can be found here.
Tip
The code below can be found in the getting-started.ipynb notebook
Step 2: Import your data¶
This example will assume the following dataset is located in your current working directory (feel free to replace it with your data): exampleData.csv
For this dataset, importing the data looks something like:
import numpy as np
data = np.genfromtxt('./exampleData.csv', delimiter=',')
frequencies = data[:,0]
Z = data[:,1] + 1j*data[:,2]
# keep only the impedance data in the first quandrant
frequencies = frequencies[np.imag(Z) < 0]
Z = Z[np.imag(Z) < 0]
Step 3: Define your impedance model¶
Next we want to define our impedance model. In order to enable a wide variety
of researchers to use the tool, impedance.py
allows you to define a
custom circuit with any combination of circuit elements.
The circuit is defined as a string (i.e. using ''
in Python), where elements in
series are separated by a dash (-
), and elements in parallel are wrapped in
a p( , )
. Each element is defined by the function (in circuit-elements.py) followed by a single digit identifier.
For example, the circuit below:

would be defined as R0-p(R1,C1)-p(R2-W1,C2)
.
Each circuit, we want to fit also needs to have an initial guess for each
of the parameters. These inital guesses are passed in as a list in order the
parameters are defined in the circuit string. For example, a good guess for this
battery data is initial_guess = [.01, .01, 100, .01, .05, 100, 1]
.
We create the circuit by importing the CustomCircuit object and passing it our circuit string and initial guesses.
from impedance.circuits import CustomCircuit
circuit = 'R0-p(R1,C1)-p(R2-W1,C2)'
initial_guess = [.01, .01, 100, .01, .05, 100, 1]
circuit = CustomCircuit(circuit, initial_guess=initial_guess)
Step 4: Fit the impedance model to data¶
Once we’ve defined our circuit, fitting it to impedance data is as easy as calling the .fit() method and passing it our experimental data!
circuit.fit(frequencies, Z)
We can access the fit parameters with circuit.parameters_
or by
printing the circuit object itself, print(circuit)
.
Step 5: Analyze/Visualize the results¶
For this dataset, the resulting fit parameters are
Parameter | Value |
\(R_0\) | 1.65e-02 |
\(R_1\) | 8.68e-03 |
\(C_1\) | 3.32e+00 |
\(R_2\) | 5.39e-03 |
\(W_{1,0}\) | 6.31e-02 |
\(W_{1,1}\) | 2.33e+02 |
\(C_2\) | 2.20e-01 |
We can get the resulting fit impedance by passing a list of frequencies to the .predict()
method.
Z_fit = circuit.predict(frequencies)
To easily visualize the fit, the plot_nyquist()
function can be handy.
import matplotlib.pyplot as plt
from impedance.plotting import plot_nyquist
fig, ax = plt.subplots()
plot_nyquist(ax, frequencies, Z, fmt='o')
plot_nyquist(ax, frequencies, Z_fit, fmt='-')
plt.legend(['Data', 'Fit'])
plt.show()

Important
🎉 Congratulations! You’re now up and running with impedance.py 🎉
Examples¶
Fitting impedance spectra¶
1. Import and initialize equivalent circuit(s)¶
To begin we will import the Randles’ circuit and a custom circuit from the impedance package. A full list of currently available circuits are available in the documentation.
[1]:
import sys
sys.path.append('../../../')
from impedance.circuits import Randles, CustomCircuit
The classes we just imported represent different equivalent circuit models. To actually use them we want to initialize a specific instance and provide an initial guess for the parameters and any other options.
E.g. for the randles circuit, one of the options is for a constant phase element (CPE) instead of an ideal capacitor.
[2]:
randles = Randles(initial_guess=[.01, .005, .1, .001, 200])
randlesCPE = Randles(initial_guess=[.01, .005, .1, .9, .001, 200], CPE=True)
Defining the custom circuit works a little differently. Here we pass a string comprised of the circuit elements grouped either in series (separated with a -
) or in parallel (using the form p(X,Y)
). Elements with multiple parameters are given in the form X/Y
[3]:
customCircuit = CustomCircuit(initial_guess=[.01, .005, .1, .005, .1, .001, 200],
circuit='R_0-p(R_1,C_1)-p(R_1,C_1)-W_1')
Each of the circuit objects we create can be printed in order to see the properties that have been defined for that circuit.
[4]:
print(randles)
-------------------------------
Circuit: Randles
Circuit string: R0-p(R1,C1)-W1
Fit: False
-------------------------------
Initial guesses:
R0 = 1.00e-02
R1 = 5.00e-03
C1 = 1.00e-01
W1_0 = 1.00e-03
W1_1 = 2.00e+02
2. Formulate data¶
Several convenience functions for importing data exist in the impedance.preprocessing module; however, here we will simply read in a .csv
file containing frequencies as well as real and imaginary impedances using the numpy package.
[5]:
import numpy as np
data = np.genfromtxt('../../../data/exampleData.csv', delimiter=',')
frequencies = data[:,0]
Z = data[:,1] + 1j*data[:,2]
# keep only the impedance data in the first quandrant
frequencies = frequencies[np.imag(Z) < 0]
Z = Z[np.imag(Z) < 0]
3. Fit the equivalent circuits to a spectrum¶
Each of the circuit classes has a .fit()
method which finds the best fitting parameters.
After fitting a circuit, the fit parameters rather that the inital guesses are shown when printing.
[6]:
randles.fit(frequencies, Z)
randlesCPE.fit(frequencies, Z)
customCircuit.fit(frequencies, Z)
print(customCircuit)
-------------------------------
Circuit: None
Circuit string: R_0-p(R_1,C_1)-p(R_1,C_1)-W_1
Fit: True
-------------------------------
Fit parameters:
R_0 = 1.65e-02 +/- 1.54e-04
R_1 = 5.31e-03 +/- 2.06e-04
C_1 = 2.32e-01 +/- 1.90e-02
R_1 = 8.77e-03 +/- 1.89e-04
C_1 = 3.28e+00 +/- 1.85e-01
W_1_0 = 6.37e-02 +/- 2.03e-03
W_1_1 = 2.37e+02 +/- 1.72e+01
4a. Predict circuit model and visualize with matplotlib¶
[7]:
import matplotlib.pyplot as plt
from impedance.plotting import plot_nyquist
f_pred = np.logspace(5,-2)
randles_fit = randles.predict(f_pred)
randlesCPE_fit = randlesCPE.predict(f_pred)
customCircuit_fit = customCircuit.predict(f_pred)
fig, ax = plt.subplots(figsize=(5,5))
plot_nyquist(ax, frequencies, Z)
plot_nyquist(ax, f_pred, randles_fit, fmt='-')
plot_nyquist(ax, f_pred, randlesCPE_fit, fmt='-')
plot_nyquist(ax, f_pred, customCircuit_fit, fmt='-')
plt.show()

4b. Or use the convenient plotting method included in the package¶
This is an experimental feature with many improvements (interactive plots, Bode plots, better confidence intervals) coming soon!!
[8]:
randles.plot(f_data=frequencies, Z_data=Z)
randlesCPE.plot(f_data=frequencies, Z_data=Z)
customCircuit.plot(f_data=frequencies, Z_data=Z)
plt.show()



[ ]:
Plotting Nyquist plots of impedance spectra¶
Plotting a basically formated Nyquist plot is as easy as 1, 2, 3…
[1]:
import matplotlib.pyplot as plt
import numpy as np
import sys
sys.path.append('../../../')
from impedance.circuits import CustomCircuit
1. Read in data¶
[2]:
data = np.genfromtxt('../../../data/exampleData.csv', delimiter=',')
frequencies = data[:,0]
Z = data[:,1] + 1j*data[:,2]
frequencies = frequencies[np.imag(Z) < 0]
Z = Z[np.imag(Z) < 0]
2. Fit a custom circuit¶
(If you want to just plot experimental data without fitting a model to it, you should check out the plotting.plot_nyquist()
function)
[3]:
circuit = CustomCircuit(initial_guess=[.01, .005, .1, .005, .1, .001, 200], circuit='R_0-p(R_1,C_1)-p(R_1,C_1)-W_1')
circuit.fit(frequencies, Z)
print(circuit)
-------------------------------
Circuit: None
Circuit string: R_0-p(R_1,C_1)-p(R_1,C_1)-W_1
Fit: True
-------------------------------
Fit parameters:
R_0 = 1.65e-02 +/- 1.54e-04
R_1 = 5.31e-03 +/- 2.06e-04
C_1 = 2.32e-01 +/- 1.90e-02
R_1 = 8.77e-03 +/- 1.89e-04
C_1 = 3.28e+00 +/- 1.85e-01
W_1_0 = 6.37e-02 +/- 2.03e-03
W_1_1 = 2.37e+02 +/- 1.72e+01
3. Plot the data and fit model with confidence bounds¶
[4]:
circuit.plot(f_data=frequencies, Z_data=Z, conf_bounds='filled')
[4]:
<matplotlib.axes._subplots.AxesSubplot at 0x115f86320>

Bonus: Easy access to all the customization of matplotlib¶
Here we plot the data, changing the size of the figure, axes label fontsize, and turning off the grid by accessing the plt.Axes() object, ax
[5]:
fig, ax = plt.subplots(figsize=(10,10))
ax = circuit.plot(ax, frequencies, Z)
ax.tick_params(axis='both', which='major', labelsize=16)
ax.grid(False)
plt.show()

Validation of EIS data¶
The Kramers-Kronig Relations¶
Electrochemical impedance spectroscopy (EIS) is built on linear systems theory which requires that the system satisfy conditions of causality, linearity, and stability. The Kramers-Kronig relations consist of a set of transformations that can be used to predict one component of the impedance from the other over the frequency limits from zero to infinity. For example, one might calculate the imaginary component of the impedance from the measured real component,
where \(Z^{\prime}(\omega)\) and \(Z^{\prime\prime}(\omega)\) are the real and imaginary components of the impedance as a function of frequency, \(\omega\). Similarly, the real part of the impedance spectrum can be calculated from the imaginary part by
The residual error between the predicted and measured impedance can then be used to determine consistency with the Kramers-Kronig relations.
Practically, however, the 0 to \(\infty\) frequency range required for integration can be difficult to measure experimentally, so several other methods have been developed to ensure Kramers-Kronig relations are met:
Measurement models¶
[1]:
import matplotlib.pyplot as plt
import numpy as np
import sys
sys.path.append('../../../')
from impedance.circuits import CustomCircuit
[2]:
# read data
data = np.genfromtxt('../../../data/exampleData.csv', delimiter=',')
f = data[:,0]
Z = data[:,1] + 1j*data[:,2]
mask = f < 1000
f = f[mask]
Z = Z[mask]
[3]:
N = 10
circuit = 'R_0'
initial_guess = [.015]
for i in range(N):
circuit += f'-p(R_{i % 9 + 1},C_{i % 9 + 1})'
initial_guess.append(.03/N)
initial_guess.append(10**(3 - 6*i/N))
meas_model = CustomCircuit(initial_guess=initial_guess, circuit=circuit)
[4]:
meas_model.fit(f, Z, method='lm')
print(meas_model.get_verbose_string())
-------------------------------
Circuit: None
Circuit string: R_0-p(R_1,C_1)-p(R_2,C_2)-p(R_3,C_3)-p(R_4,C_4)-p(R_5,C_5)-p(R_6,C_6)-p(R_7,C_7)-p(R_8,C_8)-p(R_9,C_9)-p(R_1,C_1)
Fit: True
-------------------------------
Initial guesses:
R_0 = 1.50e-02
R_1 = 3.00e-03
C_1 = 1.00e+03
R_2 = 3.00e-03
C_2 = 2.51e+02
R_3 = 3.00e-03
C_3 = 6.31e+01
R_4 = 3.00e-03
C_4 = 1.58e+01
R_5 = 3.00e-03
C_5 = 3.98e+00
R_6 = 3.00e-03
C_6 = 1.00e+00
R_7 = 3.00e-03
C_7 = 2.51e-01
R_8 = 3.00e-03
C_8 = 6.31e-02
R_9 = 3.00e-03
C_9 = 1.58e-02
R_1 = 3.00e-03
C_1 = 3.98e-03
-------------------------------
Fit parameters:
R_0 = 1.36e-02 (+/- 8.23e+05)
R_1 = 1.06e-01 (+/- 7.66e-03)
C_1 = 2.81e+03 (+/- 5.00e+01)
R_2 = 1.12e-02 (+/- 2.56e-04)
C_2 = 1.51e+03 (+/- 4.39e+01)
R_3 = 3.02e-03 (+/- 1.51e-04)
C_3 = 6.50e+02 (+/- 7.27e+01)
R_4 = 1.09e-05 (+/- 7.03e+01)
C_4 = 2.76e-02 (+/- 2.38e+05)
R_5 = 2.87e-03 (+/- 4.34e-03)
C_5 = 6.33e+01 (+/- 1.35e+02)
R_6 = 6.50e-03 (+/- 3.53e-03)
C_6 = 4.37e+00 (+/- 1.09e+00)
R_7 = 3.25e-03 (+/- 4.16e-04)
C_7 = 1.64e+00 (+/- 1.89e-01)
R_8 = 6.28e-04 (+/- 2.25e-03)
C_8 = 1.21e+02 (+/- 8.47e+02)
R_9 = 3.82e-03 (+/- 1.82e-04)
C_9 = 1.84e-01 (+/- 1.99e-02)
R_1 = 2.70e-03 (+/- 8.23e+05)
C_1 = 3.08e-07 (+/- 1.88e+02)
[5]:
from impedance.plotting import plot_nyquist
res_meas_real = (Z - meas_model.predict(f)).real/np.abs(Z)
res_meas_imag = (Z - meas_model.predict(f)).imag/np.abs(Z)
fig = plt.figure(figsize=(5,8))
gs = fig.add_gridspec(3, 1)
ax1 = fig.add_subplot(gs[:2,:])
ax2 = fig.add_subplot(gs[2,:])
# plot original data
plot_nyquist(ax1, f, Z, fmt='s')
# plot measurement model
plot_nyquist(ax1, f, meas_model.predict(f), fmt='-', scale=1e3, units='\Omega')
ax1.legend(['Data', 'Measurement model'], loc=2, fontsize=12)
# Plot residuals
ax2.plot(f, res_meas_real, '-', label=r'$\Delta_{\mathrm{Real}}$')
ax2.plot(f, res_meas_imag, '-', label=r'$\Delta_{\mathrm{Imag}}$')
ax2.set_title('Measurement Model Error', fontsize=14)
ax2.tick_params(axis='both', which='major', labelsize=12)
ax2.set_ylabel('$\Delta$ $(\%)$', fontsize=14)
ax2.set_xlabel('$f$ [Hz]', fontsize=14)
ax2.set_xscale('log')
ax2.set_ylim(-.1, .1)
ax2.legend(loc=1, fontsize=14, ncol=2)
vals = ax2.get_yticks()
ax2.set_yticklabels(['{:.0%}'.format(x) for x in vals])
plt.tight_layout()
plt.show()

The Lin-KK method¶
The lin-KK method from Schönleber et al. [1] is a quick test for checking the validity of EIS data. The validity of an impedance spectrum is analyzed by its reproducibility by a Kramers-Kronig (KK) compliant equivalent circuit. In particular, the model used in the lin-KK test is an ohmic resistor, \(R_{Ohm}\), and \(M\) RC elements.
The \(M\) time constants, \(\tau_k\), are distributed logarithmically,
and are not fit during the test (only \(R_{Ohm}\) and \(R_{k}\) are free parameters).
In order to prevent under- or over-fitting, Schönleber et al. propose using the ratio of positive resistor mass to negative resistor mass as a metric for finding the optimal number of RC elements.
The argument c
defines the cutoff value for \(\mu\). The algorithm starts at M = 3
and iterates up to max_M
until a \(\mu < c\) is reached. The default of 0.85 is simply a heuristic value based off of the experience of Schönleber et al.
If the argument c
is None
, then the automatic determination of RC elements is turned off and the solution is calculated for max_M
RC elements. This manual mode should be used with caution as under- and over-fitting should be avoided.
[1] Schönleber, M. et al. A Method for Improving the Robustness of linear Kramers-Kronig Validity Tests. Electrochimica Acta 131, 20–27 (2014) doi: 10.1016/j.electacta.2014.01.034.
[6]:
import matplotlib.pyplot as plt
import numpy as np
import sys
sys.path.append('../../../')
from impedance.validation import linKK
[7]:
# read data
data = np.genfromtxt('../../../data/exampleData.csv', delimiter=',')
f = data[:,0]
Z = data[:,1] + 1j*data[:,2]
mask = f < 1000
f = f[mask]
Z = Z[mask]
[8]:
M, mu, Z_linKK, res_real, res_imag = linKK(f, Z, max_M=100)
print('\nCompleted Lin-KK Fit\nM = {:d}\nmu = {:.2f}'.format(M, mu))
10 1.0 0.00022662281947617854
Completed Lin-KK Fit
M = 13
mu = 0.82
[9]:
from impedance.plotting import plot_nyquist
fig = plt.figure(figsize=(5,8))
gs = fig.add_gridspec(3, 1)
ax1 = fig.add_subplot(gs[:2,:])
ax2 = fig.add_subplot(gs[2,:])
# plot original data
plot_nyquist(ax1, f, Z, fmt='s')
# plot measurement model
plot_nyquist(ax1, f, Z_linKK, fmt='-', scale=1e3, units='\Omega')
ax1.legend(['Data', 'Lin-KK model'], loc=2, fontsize=12)
# Plot residuals
ax2.plot(f, res_real, '-', label=r'$\Delta_{\mathrm{Real}}$')
ax2.plot(f, res_imag, '-', label=r'$\Delta_{\mathrm{Imag}}$')
ax2.set_title('Lin-KK Model Error', fontsize=14)
ax2.tick_params(axis='both', which='major', labelsize=12)
ax2.set_ylabel('$\Delta$ $(\%)$', fontsize=14)
ax2.set_xlabel('$f$ [Hz]', fontsize=14)
ax2.set_xscale('log')
ax2.set_ylim(-.1, .1)
ax2.legend(loc=1, fontsize=14, ncol=2)
vals = ax2.get_yticks()
ax2.set_yticklabels(['{:.0%}'.format(x) for x in vals])
plt.tight_layout()
plt.show()

Circuits¶
-
class
impedance.circuits.
BaseCircuit
(initial_guess, name=None, bounds=None)[source]¶ Base class for equivalent circuit models
Methods
fit
(self, frequencies, impedance[, method, …])Fit the circuit model get_param_names
(self)Converts circuit string to names get_verbose_string
(self)Defines the pretty printing of all data in the circuit plot
(self[, ax, f_data, Z_data, …])a convenience method for plotting Nyquist plots predict
(self, frequencies[, use_initial])Predict impedance using a fit equivalent circuit model -
fit
(self, frequencies, impedance, method='lm', bounds=None)[source]¶ Fit the circuit model
Parameters: - frequencies: numpy array
Frequencies
- impedance: numpy array of dtype ‘complex128’
Impedance values to fit
- method: {‘lm’, ‘trf’, ‘dogbox’}, optional
Name of method to pass to scipy.optimize.curve_fit
- bounds: 2-tuple of array_like, optional
Lower and upper bounds on parameters. Defaults to bounds on all parameters of 0 and np.inf, except the CPE alpha which has an upper bound of 1
Returns: - self: returns an instance of self
-
plot
(self, ax=None, f_data=None, Z_data=None, conf_bounds=None, scale=1, units='Ohms')[source]¶ a convenience method for plotting Nyquist plots
Parameters: - f_data: np.array of type float
Frequencies of input data (for Bode plots)
- Z_data: np.array of type complex
Impedance data to plot
- conf_bounds: {‘error_bars’, ‘filled’, ‘filledx’, ‘filledy’}, optional
Include bootstrapped confidence bands (95%) on the predicted best fit model shown as either error bars or a filled confidence area. Confidence bands are estimated by simulating the spectra for 10000 randomly sampled parameter sets where each of the parameters is sampled from a normal distribution
Returns: - ax: matplotlib.axes
axes of the created nyquist plot
-
predict
(self, frequencies, use_initial=False)[source]¶ Predict impedance using a fit equivalent circuit model
Parameters: - frequencies: numpy array
Frequencies
- use_initial: boolean
If true and a fit was already completed, use the initial parameters instead
Returns: - impedance: numpy array of dtype ‘complex128’
Predicted impedance
-
-
class
impedance.circuits.
CustomCircuit
(circuit, **kwargs)[source]¶ Methods
fit
(self, frequencies, impedance[, method, …])Fit the circuit model get_param_names
(self)Converts circuit string to names get_verbose_string
(self)Defines the pretty printing of all data in the circuit plot
(self[, ax, f_data, Z_data, …])a convenience method for plotting Nyquist plots predict
(self, frequencies[, use_initial])Predict impedance using a fit equivalent circuit model
-
class
impedance.circuits.
Randles
(CPE=False, **kwargs)[source]¶ A Randles circuit model class
Methods
fit
(self, frequencies, impedance[, method, …])Fit the circuit model get_param_names
(self)Converts circuit string to names get_verbose_string
(self)Defines the pretty printing of all data in the circuit plot
(self[, ax, f_data, Z_data, …])a convenience method for plotting Nyquist plots predict
(self, frequencies[, use_initial])Predict impedance using a fit equivalent circuit model
Circuit Elements¶
-
impedance.circuit_elements.
A
(p, f)[source]¶ defines a semi-infinite Warburg element
Notes
\[Z = \frac{A_W}{\sqrt{ 2 \pi f}} (1-j)\]
-
impedance.circuit_elements.
C
(p, f)[source]¶ defines a capacitor
\[Z = \frac{1}{C \times j 2 \pi f}\]
-
impedance.circuit_elements.
E
(p, f)[source]¶ defines a constant phase element
Notes
\[Z = \frac{1}{Q \times (j 2 \pi f)^\alpha}\]where \(Q\) = p[0] and \(\alpha\) = p[1].
-
impedance.circuit_elements.
G
(p, f)[source]¶ defines a Gerischer Element
Notes
\[Z = \frac{1}{Y \times \sqrt{K + j 2 \pi f }}\]
-
impedance.circuit_elements.
K
(p, f)[source]¶ An RC element for use in lin-KK model
Notes
\[Z = \frac{R}{1 + j \omega \tau_k}\]
-
impedance.circuit_elements.
T
(p, f)[source]¶ A macrohomogeneous porous electrode model from Paasch et al. [1]
Notes
\[Z = A\frac{\coth{\beta}}{\beta} + B\frac{1}{\beta\sinh{\beta}}\]where
\[A = d\frac{\rho_1^2 + \rho_2^2}{\rho_1 + \rho_2} \quad B = d\frac{2 \rho_1 \rho_2}{\rho_1 + \rho_2}\]and
\[\beta = (a + j \omega b)^{1/2} \quad a = \frac{k d^2}{K} \quad b = \frac{d^2}{K}\][1] G. Paasch, K. Micka, and P. Gersdorf, Electrochimica Acta, 38, 2653–2662 (1993) doi: 10.1016/0013-4686(93)85083-B.
-
impedance.circuit_elements.
W
(p, f)[source]¶ defines a blocked boundary Finite-length Warburg Element
Notes
\[Z = \frac{R}{\sqrt{ T \times j 2 \pi f}} \coth{\sqrt{T \times j 2 \pi f }}\]where \(R\) = p[0] (Ohms) and \(T\) = p[1] (sec) = \(\frac{L^2}{D}\)
-
impedance.circuit_elements.
num_params
(n)[source]¶ decorator to store number of parameters for an element
Fitting¶
-
impedance.fitting.
buildCircuit
(circuit, frequencies, *parameters, eval_string='', index=0)[source]¶ recursive function that transforms a circuit, parameters, and frequencies into a string that can be evaluated
Parameters: - circuit: str
- parameters: list/tuple/array of floats
- frequencies: list/tuple/array of floats
Returns: - eval_string: str
Python expression for calculating the resulting fit
- index: int
Tracks parameter index through recursive calling of the function
-
impedance.fitting.
circuit_fit
(frequencies, impedances, circuit, initial_guess, method='lm', bounds=None, bootstrap=False)[source]¶ Main function for fitting an equivalent circuit to data
Parameters: - frequencies : numpy array
Frequencies
- impedances : numpy array of dtype ‘complex128’
Impedances
- circuit : string
string defining the equivalent circuit to be fit
- initial_guess : list of floats
initial guesses for the fit parameters
- method : {‘lm’, ‘trf’, ‘dogbox’}, optional
Name of method to pass to scipy.optimize.curve_fit
- bounds : 2-tuple of array_like, optional
Lower and upper bounds on parameters. Defaults to bounds on all parameters of 0 and np.inf, except the CPE alpha which has an upper bound of 1
Returns: - p_values : list of floats
best fit parameters for specified equivalent circuit
- p_errors : list of floats
one standard deviation error estimates for fit parameters
Notes
Need to do a better job of handling errors in fitting. Currently, an error of -1 is returned.
-
impedance.fitting.
computeCircuit
(circuit, frequencies, *parameters)[source]¶ evaluates a circuit string for a given set of parameters and frequencies
Parameters: - circuit : string
- frequencies : list/tuple/array of floats
- parameters : list/tuple/array of floats
Returns: - array of complex numbers
Validation¶
Interpreting EIS data fundamentally relies on the the system conforming to conditions of causality, linearity, and stability. For an example of how the adherence to the Kramers-Kronig relations, see the Validation Example Jupyter Notebook
Lin-KK method¶
Validating your data with the lin-KK model requires fitting an optimal number of RC-elements and analysis of the residual errors.
-
impedance.validation.
eval_linKK
(Rs, ts, f)[source]¶ Builds a circuit of RC elements to be used in LinKK
-
impedance.validation.
fitLinKK
(f, ts, M, Z)[source]¶ Fits the linKK model using scipy.optimize.least_squares
-
impedance.validation.
linKK
(f, Z, c=0.85, max_M=50)[source]¶ A method for implementing the Lin-KK test for validating linearity [1]
Parameters: - f: np.ndarray
measured frequencies
- Z: np.ndarray of complex numbers
measured impedances
- c: np.float
cutoff for mu
- max_M: int
the maximum number of RC elements
Returns: - mu: np.float
under- or over-fitting measure
- residuals: np.ndarray of complex numbers
the residuals of the fit at input frequencies
- Z_fit: np.ndarray of complex numbers
impedance of fit at input frequencies
Notes
The lin-KK method from Schönleber et al. [1] is a quick test for checking the validity of EIS data. The validity of an impedance spectrum is analyzed by its reproducibility by a Kramers-Kronig (KK) compliant equivalent circuit. In particular, the model used in the lin-KK test is an ohmic resistor, \(R_{Ohm}\), and \(M\) RC elements.
\[\hat Z = R_{Ohm} + \sum_{k=1}^{M} \frac{R_k}{1 + j \omega \tau_k}\]The \(M\) time constants, \(\tau_k\), are distributed logarithmically,
\[\tau_1 = \frac{1}{\omega_{max}} ; \tau_M = \frac{1}{\omega_{min}} ; \tau_k = 10^{\log{(\tau_{min}) + \frac{k-1}{M-1}\log{{( \frac{\tau_{max}}{\tau_{min}}}})}}\]and are not fit during the test (only \(R_{Ohm}\) and \(R_{k}\) are free parameters).
In order to prevent under- or over-fitting, Schönleber et al. propose using the ratio of positive resistor mass to negative resistor mass as a metric for finding the optimal number of RC elements.
\[\mu = 1 - \frac{\sum_{R_k \ge 0} |R_k|}{\sum_{R_k < 0} |R_k|}\]The argument
c
defines the cutoff value for \(\mu\). The algorithm starts atM = 3
and iterates up tomax_M
until a \(\mu < c\) is reached. The default of 0.85 is simply a heuristic value based off of the experience of Schönleber et al.If the argument
c
isNone
, then the automatic determination of RC elements is turned off and the solution is calculated formax_M
RC elements. This manual mode should be used with caution as under- and over-fitting should be avoided.[1] Schönleber, M. et al. A Method for Improving the Robustness of linear Kramers-Kronig Validity Tests. Electrochimica Acta 131, 20–27 (2014) doi: 10.1016/j.electacta.2014.01.034.
Preprocessing¶
Methods for preprocessing impedance data from instrument files
-
impedance.preprocessing.
readAutolab
(filename)[source]¶ function for reading the .csv file from Autolab
Parameters: - filename: string
Filename of .csv file to extract impedance data from
Returns: - frequencies : np.ndarray
Array of frequencies
- impedance : np.ndarray of complex numbers
Array of complex impedances
-
impedance.preprocessing.
readFile
(filename, instrument=None)[source]¶ A wrapper for reading in many common types of impedance files
Parameters: - filename: string
Filename to extract impedance data from
- instrument: string
Type of instrument file
Returns: - frequencies : np.ndarray
Array of frequencies
- impedance : np.ndarray of complex numbers
Array of complex impedances
-
impedance.preprocessing.
readGamry
(filename)[source]¶ function for reading the .DTA file from Gamry
Parameters: - filename: string
Filename of .DTA file to extract impedance data from
Returns: - frequencies : np.ndarray
Array of frequencies
- impedance : np.ndarray of complex numbers
Array of complex impedances
-
impedance.preprocessing.
readParstat
(filename)[source]¶ function for reading the .txt file from Parstat
Parameters: - filename: string
Filename of .txt file to extract impedance data from
Returns: - frequencies : np.ndarray
Array of frequencies
- impedance : np.ndarray of complex numbers
Array of complex impedances