The Inversion Class¶
Standard usage of bisip
involves instantiating an
Inversion
object. This is normally done by invoking one of the
inversion model classes described in the SIP models section below.
-
class
bisip.models.
Inversion
(filepath, nwalkers=32, nsteps=5000, headers=1, ph_units='mrad')¶ Bases:
bisip.plotlib.plotlib
,bisip.utils.utils
An abstract class to perform inversion of SIP data.
This class is the base constructor for the PolynomialDecomposition and ColeCole classes.
Parameters: - filepath (
str
) – The path to the file to perform inversion on. - nwalkers (
int
) – Number of walkers to use to explore the parameter space. Defaults to 32. - nsteps (
int
) – Number of steps to perform in the MCMC simulation. Defaults to 5000. - headers (
int
) – The number of header lines in the file. Defaults to 1. - ph_units (
str
) – The units of the phase shift measurements. Choices: ‘mrad’, ‘rad’, ‘deg’. Defaults to ‘mrad’.
-
data
¶ The input data dictionary.
Type: dict
-
fit
(p0=None, pool=None, moves=None)¶ Samples the posterior distribution to fit the model to the data.
Parameters: - p0 (
ndarray
) – Starting parameter values. Should be a 2D array with shape (nwalkers, ndim). If None, random values will be uniformly drawn from the parameter bounds. Defaults to None. - pool (
pool
, optional) – A pool object from the Python multiprocessing library. See https://emcee.readthedocs.io/en/stable/tutorials/parallel/. Defaults to None. - moves (
moves
, optional) – A emcee Moves class (see https://emcee.readthedocs.io/en/stable/user/moves/). If None, the emcee algorithm StretchMove is used. Defaults to None.
- p0 (
-
fitted
¶ Whether the model has been fitted or not.
Type: bool
-
get_chain
(**kwargs)¶ Gets the MCMC chains from a fitted model.
Keyword Arguments: - discard (
int
) – Number of steps to discard (burn-in period). - thin (
int
) – Thinning factor. - flat (
bool
) – Whether or not to flatten the walkers. If flat is False, the output chain will have shape (nsteps, nwalkers, ndim). If flat is True, the output chain will have shape (nsteps*nwalkers, ndim).
Returns: The MCMC chain(s).
Return type: ndarray
- discard (
-
p0
¶ Starting parameter values. Should be a 2D array with shape (nwalkers, ndim).
Type: ndarray
-
param_bounds
¶ Ordered bounds of the parameters.
Type: list
offloat
-
param_names
¶ Ordered names of the parameters.
Type: list
ofstr
-
params
¶ Parameter names and their bounds.
Type: dict
-
sampler
¶ A emcee sampler object (see https://emcee.readthedocs.io/en/stable/user/sampler/).
Type: EnsembleSampler
- filepath (
SIP models¶
Pelton Cole-Cole¶
-
class
bisip.models.
PeltonColeCole
(*args, n_modes=1, **kwargs)¶ Bases:
bisip.models.Inversion
A generalized ColeCole inversion scheme for SIP data.
Parameters: - *args – Arguments passed to the Inversion class.
- n_modes (
int
) – The number of ColeCole modes to use for the inversion. Defaults to 1. - **kwargs – Additional keyword arguments passed to the Inversion class.
-
forward
(theta, w)¶ Returns a ColeCole impedance.
Parameters: - theta (
ndarray
) – Ordered array of R0, m_{1}, …, m_{n_modes}, log_tau_{1}, …, log_tau_{n_modes}, c_{1}, …, c_{n_modes}. See https://doi.org/10.1016/j.cageo.2017.05.001. - w (
ndarray
) – Array of angular frequencies to compute the impedance for (w = 2*pi*f).
- theta (
Polynomial decomposition¶
-
class
bisip.models.
PolynomialDecomposition
(*args, poly_deg=5, c_exp=1.0, **kwargs)¶ Bases:
bisip.models.Inversion
A polynomial decomposition inversion scheme for SIP data.
Parameters: - *args – Arguments passed to the Inversion class.
- poly_deg (
int
) – The polynomial degree to use for the decomposition. Defaults to 5. - c_exp (
float
) – The c-exponent to use for the decomposition scheme. 0.5 -> Warburg, 1.0 -> Debye. Defaults to 1.0. - **kwargs – Additional keyword arguments passed to the Inversion class.
-
forward
(theta, w)¶ Returns a Polynomial Decomposition impedance.
Parameters: - theta (
ndarray
) – Ordered array of R0, a_{poly_deg}, a_{poly_deg-1}, …, a_{0}. See https://doi.org/10.1016/j.cageo.2017.05.001. - w (
ndarray
) – Array of angular frequencies to compute the impedance for (w = 2*pi*f).
- theta (
Dias (2000)¶
-
class
bisip.models.
Dias2000
(*args, **kwargs)¶ Bases:
bisip.models.Inversion
A generalized ColeCole inversion scheme for SIP data.
Parameters: - *args – Arguments passed to the Inversion class.
- **kwargs – Additional keyword arguments passed to the Inversion class.
-
forward
(theta, w)¶ Returns a Dias (2000) impedance.
Parameters: - theta (
ndarray
) – Ordered array of R0, m, log_tau, eta, delta. See https://doi.org/10.1016/j.cageo.2017.05.001. - w (
ndarray
) – Array of angular frequencies to compute the impedance for (w = 2*pi*f).
- theta (
Shin (2015)¶
-
class
bisip.models.
Shin2015
(*args, **kwargs)¶ Bases:
bisip.models.Inversion
A Shin (2015) inversion scheme for SIP data.
Parameters: - *args – Arguments passed to the Inversion class.
- **kwargs – Additional keyword arguments passed to the Inversion class.
Warning
The Shin model implementation is yielding unexpected results and needs to be reviewed.
-
forward
(theta, w)¶ Returns a Shin (2015) impedance.
Parameters: - theta (
ndarray
) – Ordered array of R1, R2, log_Q1, log_Q2, n1, n2. - w (
ndarray
) – Array of angular frequencies to compute the impedance for (w = 2*pi*f).
- theta (
Plotting methods¶
These functions may be called as methods of the Inversion
class
after fitting the model to a dataset.
-
class
bisip.plotlib.
plotlib
¶ -
plot_corner
(chain=None, **kwargs)¶ Plots the corner plot of the MCMC simulation.
Parameters: - chain (
ndarray
) – A numpy array containing the MCMC chain to plot. Should have a shape (nsteps, nwalkers, ndim) or (nsteps*nwalkers, ndim). If None, the full, unflattened chain will be used and all walkers will be plotted. Defaults to None. - **kwargs – Additional keyword arguments for the get_chain function (see below). Use these arguments only if not explicitly passing the chain argument.
Keyword Arguments: - discard (
int
) – The number of steps to discard. - thin (
int
) – The thinning factor (keep every thin step). - flat (
bool
) – Whether to flatten the walkers into a single chain or not.
Returns: A matplotlib figure.
Return type: Figure
- chain (
-
plot_data
(feature='phase', **kwargs)¶ Plots the input data.
Parameters: - feature (
str
) – Which data feature to plot. Choices are ‘phase’, ‘amplitude’, ‘real’, ‘imaginary’. Defaults to ‘phase’. - **kwargs – Additional keyword arguments passed to the matplotlib errorbar function.
Returns: A matplotlib figure.
Return type: Figure
- feature (
-
plot_fit
(chain=None, p=[2.5, 50, 97.5], **kwargs)¶ Plots the input data, best fit and confidence interval of a model. Shows the real and imaginary parts.
Parameters: - chain (
ndarray
) – A numpy array containing the MCMC chain to plot. Should have a shape (nsteps, nwalkers, ndim) or (nsteps*nwalkers, ndim). If None, the full, unflattened chain will be used and all walkers will be plotted. Defaults to None. - p (
list
ofint
) – Percentile values for lower confidence interval, best fit curve, and upper confidence interval, in that order. Defaults to [2.5, 50, 97.5] for the median and 95% HPD. - **kwargs – Additional keyword arguments for the get_chain function (see below). Use these arguments only if not explicitly passing the chain argument.
Keyword Arguments: - discard (
int
) – The number of steps to discard. - thin (
int
) – The thinning factor (keep every thin step). - flat (
bool
) – Whether to flatten the walkers into a single chain or not.
Returns: A matplotlib figure.
Return type: Figure
- chain (
-
plot_fit_pa
(chain=None, p=[2.5, 50, 97.5], **kwargs)¶ Plots the input data, best fit and confidence interval of a model. Shows the amplitude and phase spectra.
Parameters: - chain (
ndarray
) – A numpy array containing the MCMC chain to plot. Should have a shape (nsteps, nwalkers, ndim) or (nsteps*nwalkers, ndim). If None, the full, unflattened chain will be used and all walkers will be plotted. Defaults to None. - p (
list
ofint
) – Percentile values for lower confidence interval, best fit curve, and upper confidence interval, in that order. Defaults to [2.5, 50, 97.5] for the median and 95% HPD. - **kwargs – Additional keyword arguments for the get_chain function (see below). Use these arguments only if not explicitly passing the chain argument.
Keyword Arguments: - discard (
int
) – The number of steps to discard. - thin (
int
) – The thinning factor (keep every thin step). - flat (
bool
) – Whether to flatten the walkers into a single chain or not.
Returns: A matplotlib figure.
Return type: Figure
- chain (
-
plot_histograms
(chain=None, bins=25, **kwargs)¶ Plots histograms of the MCMC simulation chains.
Parameters: - chain (
ndarray
) – A numpy array containing the MCMC chain to plot. Should have a shape (nsteps, nwalkers, ndim) or (nsteps*nwalkers, ndim). If None, the full, unflattened chain will be used and all walkers will be plotted. Defaults to None. - bins (
int
) – The number of bins to use in the histograms. - **kwargs – Additional keyword arguments for the get_chain function (see below). Use these arguments only if not explicitly passing the chain argument.
Keyword Arguments: - discard (
int
) – The number of steps to discard. - thin (
int
) – The thinning factor (keep every thin step). - flat (
bool
) – Whether to flatten the walkers into a single chain or not.
Returns: A matplotlib figure.
Return type: Figure
- chain (
-
plot_traces
(chain=None, **kwargs)¶ Plots the traces of the MCMC simulation.
Parameters: - chain (
ndarray
) – A numpy array containing the MCMC chain to plot. Should have a shape (nsteps, nwalkers, ndim) or (nsteps*nwalkers, ndim). If None, the full, unflattened chain will be used and all walkers will be plotted. Defaults to None. - **kwargs – Additional keyword arguments for the get_chain function (see below). Use these arguments only if not explicitly passing the chain argument.
Keyword Arguments: - discard (
int
) – The number of steps to discard. - thin (
int
) – The thinning factor (keep every thin step). - flat (
bool
) – Whether to flatten the walkers into a single chain or not.
Returns: A matplotlib figure.
Return type: Figure
- chain (
-
Utility methods¶
These utility functions may be called as methods of the Inversion
class.
-
class
bisip.utils.
utils
¶ -
get_model_percentile
(p=[2.5, 50, 97.5], chain=None, **kwargs)¶ Gets percentiles of the model values for a MCMC chain.
Parameters: - p (
float
orlist
offloat
) – percentiles values to compute. Defaults to [2.5, 50, 97.5]. - chain (
ndarray
) – A numpy array containing the MCMC chain to plot. Should be a 2D array (nsteps, ndim). If None and no kwargs are passed to discard iterations, will raise a warning and the full chain will be used. Defaults to None.
Keyword Arguments: **kwargs – See kwargs of the get_chain method.
- p (
-
get_param_mean
(chain=None, **kwargs)¶ Gets the mean of the model parameters for a MCMC chain.
Parameters: chain ( ndarray
) – A numpy array containing the MCMC chain to plot. Should be a 2D array (nsteps, ndim). If None and no kwargs are passed to discard iterations, will raise a warning and the full chain will be used. Defaults to None.Keyword Arguments: **kwargs – See kwargs of the get_chain method.
-
get_param_percentile
(p=[2.5, 50, 97.5], chain=None, **kwargs)¶ Gets percentiles of the parameter values for a MCMC chain.
Parameters: - p (
float
orlist
offloat
) – percentiles values to compute. Defaults to [2.5, 50, 97.5]. - chain (
ndarray
) – A numpy array containing the MCMC chain to plot. Should be a 2D array (nsteps, ndim). If None and no kwargs are passed to discard iterations, will raise a warning and the full chain will be used. Defaults to None.
Keyword Arguments: **kwargs – See kwargs of the get_chain method.
- p (
-
get_param_std
(chain=None, **kwargs)¶ Gets the standard deviation of the model parameters.
Parameters: chain ( ndarray
) – A numpy array containing the MCMC chain to plot. Should be a 2D array (nsteps, ndim). If None and no kwargs are passed to discard iterations, will raise a warning and the full chain will be used. Defaults to None.Keyword Arguments: **kwargs – See kwargs of the get_chain method.
-
load_data
(filename, headers=1, ph_units='mrad')¶ Imports a data file and prepares it for inversion.
Parameters: - filepath (
str
) – The path to the data file. - headers (
int
) – The number of header lines in the file. Defaults to 1. - ph_units (
str
) – The units of the phase shift measurements. Choices: ‘mrad’, ‘rad’, ‘deg’. Defaults to ‘mrad’.
- filepath (
-
print_latex_parameters
(names, values, uncertainties, decimals=3)¶ Prints Pelton parameters and their uncertainties with LaTeX.
Requires the code to be run in interactive mode either from a Jupyter notebook or IPython console to use display and Math functions.
Parameters: - names (
list
) – The list of parameter names. Obtain it fromInversion.param_names
. - values (
list
) – Ordered list of mean parameter values. Obtain it fromInversion.get_param_mean()
. - uncertainties (
list
) – Ordered list of parameter uncertainties. Obtain it fromInversion.get_param_std()
. - decimals (
int
) – The number of decimals to display.
- names (
-