API
cfg
Package-wide settings and constants.
- pysoplot.cfg.lam238 = 0.000155125479614
238U decay constant [Ma^-1]. Default [JAFFEY1971].
- pysoplot.cfg.lam235 = 0.000984849860843
235U decay constant [Ma^-1]. Default [JAFFEY1971].
- pysoplot.cfg.lam234 = 2.822030700105632
234U decay constant [Ma^-1]. Default [CHENG2013].
- pysoplot.cfg.lam231 = 21.15511004303205
231Pa decay constant [Ma^-1]. Default [ROBERT1969].
- pysoplot.cfg.lam230 = 9.170554357535263
230Th decay constant [Ma^-1]. Default [CHENG2013].
- pysoplot.cfg.lam226 = 433.2169878499658
226Ra decay constant [Ma^-1]
- pysoplot.cfg.lam232 = 4.947517348750502e-05
232Th decay constant [Ma^-1] (default [HOLDEN1990])
- pysoplot.cfg.s238 = 8.332053601458737e-08
238U decay constant 1 sigma uncertainty [Ma^-1]. Default [JAFFEY1971].
- pysoplot.cfg.s235 = 6.716698160081028e-07
235U decay constant 1 sigma uncertainty [Ma^-1]. Default [JAFFEY1971].
- pysoplot.cfg.s234 = 0.001493624261109568
234U decay constant 1 sigma uncertainty [Ma^-1]. Default [CHENG2013].
- pysoplot.cfg.s231 = 0.07102280191465059
231Pa decay constant 1 sigma uncertainty [Ma^-1]. Default [ROBERT1969].
- pysoplot.cfg.s230 = 0.006673111897550267
230Th decay constant 1 sigma uncertainty [Ma^-1]. Default [CHENG2013].
- pysoplot.cfg.s226 = 1.8953243218436007
226Ra decay constant 1 sigma uncertainty [Ma^-1]
- pysoplot.cfg.s232 = 3.531422306388949e-07
232Th decay constant 1 sigma uncertainty [Ma^-1]. Default [HOLDEN1990].
- pysoplot.cfg.U = 137.818
Modern natural 238U/235U ratio. Default [HIESS2012].
- pysoplot.cfg.sU = 0.0225
Modern natural 238U/235U ratio 1 sigma uncertainty. Default [HIESS2012].
- pysoplot.cfg.rng
Random number generator used across all modules. Allows for reproducible Monte Carlo simulations results.
- Type
numpy.Generator
- pysoplot.cfg.h = 1.4
Spine h value (see [Powell2020]).
- Type
float
- pysoplot.cfg.mswd_ci_thresholds = [0.85, 0.95]
MSWD one-sided confidence interval thresholds for classical regression fits. First element is the model 1 –> 1x threshold, second element is the model 1x –> 2/3 threshold.
- pysoplot.cfg.mswd_wav_ci_thresholds = [0.85, 0.95]
MSWD one-sided confidence interval thresholds for classical wtd. averages. First element is analytical –> analytical + observed scatter threshold. Second element is not yet used. Equivalent Isoplot default for first element is 0.70.
- pysoplot.cfg.axis_labels_kw
Axis labels key-word arguments.
- pysoplot.cfg.conc_age_ellipse_kw
Concordia age ellipse marker key-word arguments.
- pysoplot.cfg.conc_env_kw
Concordia uncertainty envelope fill key-word arguments.
- pysoplot.cfg.conc_env_line_kw
Concordia envelope line key-word arguments.
- pysoplot.cfg.conc_intercept_markers_kw
Concordia intercept markers key-word arguments.
- pysoplot.cfg.conc_line_kw
Concordia line key-word arguments.
- pysoplot.cfg.conc_markers_kw
Concordia age markers key-word arguments.
- pysoplot.cfg.conc_text_kw
Concordia marker labels key-word arguments. Caution: Be careful changing the annotation_clip and clip_on settings.
- pysoplot.cfg.dp_ellipse_kw
Data point confidence ellipse key-word arguments.
- pysoplot.cfg.dp_label_kw
Data point confidence ellipse key-word arguments.
- pysoplot.cfg.fig_kw
Figure key-word arguments.
- pysoplot.cfg.gridlines_kw
Figure gridline key-word arguments.
- pysoplot.cfg.hist_bars_kw
Histogram bar key-word arguments.
- pysoplot.cfg.hist_fig_kw
Histogram figure key-word arguments.
- pysoplot.cfg.major_ticks_kw
Major axis tick key-word arguments.
- pysoplot.cfg.minor_ticks_kw
Minor axis tick key-word arguments.
- pysoplot.cfg.pb76_line_kw
Common 207Pb/206Pb projection line key-word arguments.
- pysoplot.cfg.renv_kw
Regression envelope key-word arguments.
- pysoplot.cfg.renv_line_kw
Regression envelope line key-word arguments.
- pysoplot.cfg.rline_kw
Regression line key-word arguments.
- pysoplot.cfg.scatter_markers_kw
Scatter plot marker key-word arguments.
- pysoplot.cfg.spine_kw
Figure axis spine key-word arguments. Spines are the axis lines.
- pysoplot.cfg.subplot_kw
subplot_kw passed to matplotlib.pyplot.subplots().
- pysoplot.cfg.wav_env_kw
Weighted average envelope key-word arguments.
- pysoplot.cfg.wav_fig_kw
Weighted average figure key-word arguments.
- pysoplot.cfg.wav_line_kw
Weighted average line key-word arguments.
- pysoplot.cfg.wav_markers_kw
Weighted average data point marker key-word arguments.
concordia
Concordia plotting routines for (dis)equilibrium U-Pb datasets.
- pysoplot.concordia.plot_concordia(ax, diagram, plot_markers=True, env=False, age_ellipses=False, marker_max=None, marker_ages=(), auto_markers=True, remove_overlaps=True, age_prefix='Ma')
Plot equilibrium U-Pb concordia curve on concordia diagram.
- Parameters
ax (matplotlib.pyplot.Axes) – Axes object to plot concordia curve in.
diagram ({'tw', 'wc'}) – Concordia diagram type.
marker_max (float, optional) – User specified age marker max (Ma).
- Raises
UserWarning – if concordia lies entirely outside the axis limits.:
- pysoplot.concordia.plot_diseq_concordia(ax, A, init, diagram, sA=None, age_ellipses=False, plot_markers=True, marker_max=None, env=False, env_method='mc', marker_ages=(), auto_markers=True, spaghetti=False, pA=None, remove_overlaps=True, env_trials=1000, negative_ratios=True, age_prefix='Ma')
Plot disequilibrium U-Pb concordia curve on concordia diagram.
- Parameters
ax (matplotlib.pyplot.Axes) – Axes object to plot concordia curve in.
A (array-like) – one-dimensional array of activity ratio values arranged as follows - [234U/238U], [230Th/238U], [226Ra/238U], [231Pa/235U]
init (array-like) – two-element list of boolean values, the first is True if [234U/238U] is an initial value and False if a present-day value, the second is True if [230Th/238U] is an initial value and False if a present-day value
sA (array-like, optional) – one-dimensional array of activity ratio value uncertainties given as 1 sigma absolute and arranged in the same order as A
diagram ({'tw', 'wc'}) – Concordia diagram type.
marker_max (float, optional) – User specified age marker max (Ma).
spaghetti (bool) – plot each simulated line using arbitrary colours
- Raises
UserWarning – if concordia lies entirely outside the axis limits.:
dqpb
Functions and routines for computing disequilibrium U-Pb ages.
- pysoplot.dqpb.concint_age(fit, A, sA, init, t0, diagram='tw', dc_errors=False, trials=50000, u_errors=False, negative_ratios=True, negative_ages=True, intercept_plot=True, hist=(False, False), conc_kw=None, intercept_plot_kw=None, A48i_lim=(0.0, 20.0), A08i_lim=(0.0, 10.0), age_lim=(0.0, 20.0), uncert='mc')
Compute disequilibrium U-Pb concordia intercept age and age uncertainties using Monte Carlo simulation. Optionally, produce a plot of the concordia intercept.
Notes
Only accepts data in 2-D Tera-Wasserburg form at present.
- Parameters
fit (dict) – Linear regression fit parameters.
A (array-like) – One-dimensional array of activity ratio values arranged as follows - [234U/238U], [230Th/238U], [226Ra/238U], [231Pa/235U].
sA (array-like) – One-dimensional array of activity ratio value uncertainties given as 1 sigma absolute and arranged in the same order as A.
init (array-like) – Two-element list of boolean values, the first is True if [234U/238U] is an initial value and False if a present-day value, the second is True if [230Th/238U] is an initial value and False if a present-day value.
t0 (float) – Initial guess for numerical age solving (the equilibrium age is typically good enough).
uncert ({'mc', 'none'}) – Method of computing age uncertainties.
diagram ({'tw'}) – Concordia diagram type. Currently only accepts ‘tw’, i.e. Tera-Waserburg.
- pysoplot.dqpb.forced_concordance(fit57, fit86, A, sA, t0=1.0, norm_isotope='204Pb', negative_ratios=True, negative_ages=True, hist=(False, False), trials=50000)
Compute “forced concordance” [234U/238U] value following the approach of Engel et al. (2019).
- Parameters
fit57 (dict) – 207Pb isochron linear regression fit
fit86 (dict) – 206Pb isochron linear regression fit
A (array-like) – one-dimensional array of activity ratio values arranged as follows - np.nan, [230Th/238U], [226Ra/238U], [231Pa/235U]
sA (array-like) – one-dimensional array of activity ratio value uncertainties given as 1 sigma absolute and arranged in the same order as A
References
Engel, J., Woodhead, J., Hellstrom, J., Maas, R., Drysdale, R., Ford, D., 2019. Corrections for initial isotopic disequilibrium in the speleothem U-Pb dating method. Quaternary Geochronology 54, 101009. https://doi.org/10.1016/j.quageo.2019.101009
- pysoplot.dqpb.isochron_age(fit, A, sA, t0, init=(True, True), age_type='iso-206Pb', norm_isotope='204Pb', hist=(False, False), trials=50000, dc_errors=False, negative_ratios=True, negative_ages=True)
Compute disequilibrium 238U-206Pb or 235U-207Pb isochron age and Monte Carlo age uncertainties.
- Parameters
fit (dict) – linear regression fit parameters
A (array-like) – one-dimensional array of activity ratio values arranged as follows - [234U/238U], [230Th/238U], [226Ra/238U], [231Pa/235U]
sA (array-like) – one-dimensional array of activity ratio value uncertainties given as 1 sigma absolute and arranged in the same order as A
init (array-like) – two-element list of boolean values, the first is True if [234U/238U] is an initial value and False if a present-day value, the second is True if [230Th/238U] is an initial value and False if a present-day value
t0 (float) – initial guess for numerical age solving (the equilibrium age is typically good enough)
age_type ({'iso-206Pb', 'iso-207Pb'}) – Isochron age type.
- pysoplot.dqpb.pbu_age(x, Vx, t0, DThU=None, DThU_1s=None, DPaU=None, DPaU_1s=None, alpha=None, alpha_1s=None, age_type='206Pb*', uncert='analytical', rand=False, wav=False, wav_opts=None, mc_opts=None)
Compute disequilibrium radiogenic 206Pb*/238U, 207Pb*/235U ages, or 207Pb-corrected age assuming a constant ratio of minerl-melt partition coefficients, and optionally compute a weighted average age.
- Parameters
x (ndarray (1-D, 2-D)) – Either a 1 x n array of measured 206Pb*/238U for each aliquot (for 206Pb* age), or 2 x n array of measured 238U/206Pb and 207Pb/206Pb (for 207Pb-corrected age).
Vx (ndarray (2-D)) – Covariance matrix of uncertainties on measured isotope ratios. This should be an n x n array for 206Pb* and 207Pb* ages, or a 2n x 2n array for 207Pb-corrected ages.
DThU (float, optional) – Ratio of mineral-melt distribution coefficients, DTh/DU.
DThU_1s (float, optional) – Uncertainty on ratio of mineral-melt distribution coefficients, DTh/DU.
DPaU (float, optional) – Ratio of mineral-melt distribution coefficients, DPu/DU.
DPaU_1s (float, optional) – Uncertainty on ratio of mineral-melt distribution coefficients, DPu/DU.
alpha (float, optional) – Common 207Pb/206Pb ratio.
alpha_1s (float, optional) – Uncertainty on common 207Pb/206Pb ratio.
t0 (ndarray or float, optional) – Age guess(es).
uncert ({'mc', 'analytical'}, optional) – Method of propagating age uncertainties.
rand (bool) – If true, compute random only uncertainties as well as total uncertainties (e.g., for plotting error bars).
mc_opts (dict) – Monte Carlo simulation options.
wav_opts (dict) – Weighted average calculation and plotting options.
- Returns
results – Age and uncertainty results.
- Return type
dict
- pysoplot.dqpb.pbu_iterative_age(x, Vx, ThU_melt, ThU_melt_1s, t0, Pb208_206=None, V_Pb208_206=None, Th232_U238=None, V_Th232_U238=None, DPaU=None, DPaU_1s=None, alpha=None, alpha_1s=None, age_type='206Pb*', uncert='analytical', rand=False, wav=False, wav_opts=None, mc_opts=None)
Compute disequilibrium 206Pb*/U238 or 207Pb-corrected ages iteratively along with age uncertainties. Either initial Th/U_min can be inferred from measured 232Th/238U, or from radiogenic 208Pb/206Pb and age.
- Parameters
x (ndarray (1-D, 2-D)) – Either a 1 x n array of measured 206Pb*/238U for each aliquot (for 206Pb* age), or 2 x n array of measured 238U/206Pb and 207Pb/206Pb (for 207Pb-corrected age).
Vx (ndarray (2-D)) – Covariance matrix of uncertainties on measured isotope ratios. This should be an n x n array for 206Pb* and 207Pb* ages, or a 2n x 2n array for 207Pb-corrected ages.
ThU_melt (float) – Th/U ratio of the melt.
ThU_melt_1s (float) – Uncertainty on Th/U ratio of the melt (1 sigma).
Th232_U238 (ndarray (1-D), optional) – Measured 232Th/238U for each aliquot.
V_Th232_U238 (ndarray (2-D), optional) – Covariance matrix of uncertainties on measured 232Th/238U.
Pb208_206 (ndarray (1-D)) – Measured radiogenic 208Pb/206Pb for each aliquot.
V_Pb208_206 (ndarray (2-D):) – Covariance matrix of uncertainties on radiogenic 208Pb/206Pb values.
DPaU (float, optional) – Ratio of mineral-melt distribution coefficients, DPu/DU.
DPaU_1s (float, optional) – Uncertainty on ratio of mineral-melt distribution coefficients, DPu/DU.
alpha (float, optional) – Common 207Pb/206Pb ratio.
alpha_1s (float, optional) – Uncertainty on common 207Pb/206Pb ratio.
t0 (ndarray or float, optional) – Age guess(es).
uncert ({'mc', 'analytical'}, optional) – Method of propagating age uncertainties.
rand (bool) – If true, compute random only uncertainties as well as total uncertainties (e.g., for plotting error bars).
mc_opts (dict) – Monte Carlo simulation options.
wav_opts (dict) – Weighted average calculation and plotting options.
- Returns
results – Age and uncertainty results.
- Return type
dict
exceptions
Custom exceptions.
- exception pysoplot.exceptions.ConvergenceError
Exception raised when a numerical routine fails to meet convergence criterion.
isochron
“Classical” isochron ages.
- pysoplot.isochron.age(b, sb=None, age_type='iso-206Pb', dc_errors=False)
Calculate classical isochron age.
- Parameters
b (float) – Linear regression slope
sb (float, optional) – Uncertaintiy in regression slope (\(1\sigma\))
age_type ({'iso-Pb6U8', 'iso-Pb7U5'}) – Isochron age type.
- pysoplot.isochron.mc_uncert(fit, age_type='iso-206Pb', dc_errors=False, norm_isotope='204Pb', trials=50000, hist=False)
Compute classical isochron age uncertainties using Monte Carlo approach.
Compute Monte Carlo age uncertainties for equilibrium concordia intercept age.
- Parameters
fit (dict) – Linear regression fit parameters.
age_type ({'iso-Pb6U8', 'iso-Pb7U5'}) – Isochron age type.
trials (int) – Number of Monte Carlo trials.
plotting
Plotting routines and functions.
- pysoplot.plotting.apply_plot_settings(fig, plot_type='isochron', diagram=None, xlim=(None, None), ylim=(None, None), axis_labels=(None, None), norm_isotope='204Pb')
Set axis labels and limits. Apply label, tick, grid formatting settings as defined in pysoplot.cfg.
- Parameters
fig (matplotlib.pyplot.Figure) – Figure to apply settings to.
plot_type ({'isochron', 'intercept', 'wav', 'hist'}, optional) – Plot type.
- pysoplot.plotting.plot_dp(x, sx, y, sy, r_xy, labels=None, p=0.95, reset_axis_limits=True)
Plot 2-dimensional data points with assigned uncertainties as confidence ellipses.
- Parameters
x (np.ndarray) – x values (as 1-dimensional array)
sx (np.ndarray) – uncertainty on x at \(1\sigma\) abs.
y (np.ndarray) – y values
sy (np.ndaray) – uncertainty on y at \(1\sigma\) abs.
r_xy (np.ndarray) – x-y correlation coefficient
Notes
Assumes a large sample size was used to estimate data point uncertainties and covariances.
- pysoplot.plotting.plot_rfit(ax, fit)
Plot regression line and uncertainty envelope about regression line indicating 95% confidence limits on fit.
- Parameters
ax (matplotlib.pyplot.Axes) – Axes to plot fit to.
fit (dict) – Linear regression fit results.
- pysoplot.plotting.renv(ax, fit)
Plot uncertainty envelope about regression line indicating 95% confidence limits on fit. For classical regression fits, uses the approach of Ludwig (1980), for non-classical fits uses Monte Carlo simulation.
- Parameters
ax (matplotlib.pyplot.Axes) – Axes to plot to.
fit (dict) – Linear regression fit results.
References
Ludwig, K.R., 1980. Calculation of uncertainties of U-Pb isotope data. Earth and Planetary Science Letters 212–202.
- pysoplot.plotting.rline(ax, theta)
Plot regression line.
- Parameters
ax (matplotlib.pyplot.Axes) – Axes to plot to.
theta (array-like) – Array with y-intercept as first element and slope as second.
- pysoplot.plotting.wav_plot(x, xpm, xb, xbpm, rand_pm=None, sorted=False, ylim=(None, None), x_multiplier=1.0, dp_labels=None)
Plot weighted average as line and uncertainty band and data points as uncertainty bars.
- Parameters
x (np.ndarray) – Values to be averaged as 1-dimensional array.
xpm (np.ndarray) – Symmetrical +/- errors on x.
xb (float) – Average value (x-bar).
xbpm (float) – Symmetrical +/- error on xb.
rand_pm (np.ndarray, optional) – Random symmetrical +/- errors on x.
x_multiplier – Use for unit conversion (e.g. set to 1000 for Ma to ka conversion)
dp_labels (array-like, optional) – Data point labels.
regression
Linear regression algorithms for 2-dimensional datasets.
- pysoplot.regression.classical_fit(x, sx, y, sy, r_xy, model='ca', plot=False, diagram=None, isochron=True, xlim=(None, None), ylim=(None, None), axis_labels=(None, None), norm_isotope=None, dp_labels=None)
Fit a classical linear regression line to 2-dimensional data and optionally create plot of the fit. If model is set to ‘ca’, then routine will emulate default the protocols of Isoplot Ludwig (2012).
- Parameters
x (np.ndarray) – x values (as 1-dimensional array)
sx (np.ndarray) – analytical uncertainty on x at \(1\sigma\) abs.
y (np.ndarray) – y values
sy (np.ndaray) – analytical uncertainty on y at \(1\sigma\) abs.
r_xy (np.ndarray) – x-y correlation coefficient
isochron (bool, optional) – If true and model is ‘ca’ then fits model 3 for excess scatter data sets instead of model 2.
model ({'ca', 'c1', 'c2', 'c3'}) – Regression model to fit to data.
Notes
model should be one of
‘ca’ fit ‘best’ model depending on MSWD of york fit; emulates Isoplot behaviour
‘c1’ standard York fit with analytical errors
‘c2’ model of McSaveney in Faure (1977), equivalent to Isoplot model 2
‘c3’ equivalent to Isoplot model 3
References
Faure, G., 1977. Appendix 1: Fitting of isochrons for dating by the Rb-Sr method, in: Principles of Isotope Geology. John Wiley and Sons, pp. 1–17.
Ludwig, K.R., 2012. Isoplot/Ex Version 3.75: A Geochronological Toolkit for Microsoft Excel, Special Publication 4. Berkeley Geochronology Center.
- pysoplot.regression.lad(x, y)
Fit least absolute deviation model of Sadovski (1974).
Notes
Based on code by Roger Powell.
References
Sadovski, A.N., 1974. Algorithm AS 74: L1-norm Fit of a Straight Line. Journal of the Royal Statistical Society. Series C (Applied Statistics) 23, 244–248. https://doi.org/10.2307/2347013
- pysoplot.regression.lsq(x, y)
Fit ordinary least-squares model
Notes
Based on code by Roger Powell.
- pysoplot.regression.model_3(x, sx, y, sy, r_xy, sy_excess0, theta0, mswd_tol=0.0001, itmax=100, york_itmax=50, york_atol=1e-08, york_rtol=1e-08)
Fit a model 3 regression line that is equivalent to Isoplot model 3 Ludwig (2012). Weights each point according to analytical errors, plus an extra component of Gaussian distributed scatter in y (applied equally to each data point). For each iteration this excess error in y is inflated and a new York fit produced, until MSWD converges to 1.
- Parameters
x (np.ndarray) – x values (as 1-dimensional array)
sx (np.ndarray) – analytical uncertainty on x at \(1\sigma\) abs.
y (np.ndarray) – y values
sy (np.ndaray) – analytical uncertainty on y at \(1\sigma\) abs.
r_xy (np.ndarray) – x-y correlation coefficient
Notes
In Isoplot, this is only used for classical isochron ages
When applying this algorithm,there should be a good justification for the assumption of Guassian distributed excess scatter.
References
Ludwig, K.R., 2012. Isoplot/Ex Version 3.75: A Geochronological Toolkit for Microsoft Excel, Special Publication 4. Berkeley Geochronology Center.
- pysoplot.regression.robust_fit(x, sx, y, sy, r_xy, model='ra', plot=False, diagram=None, xlim=(None, None), ylim=(None, None), axis_labels=(None, None), norm_isotope='204Pb', dp_labels=None)
Fit a robust linear regression to 2-dimensional data and optionally create a plot of the fit.
- Parameters
x (np.ndarray) – x values (as 1-dimensional array)
sx (np.ndarray) – analytical uncertainty on x at \(1\sigma\) abs.
y (np.ndarray) – y values
sy (np.ndaray) – analytical uncertainty on y at \(1\sigma\) abs.
r_xy (np.ndarray) – x-y correlation coefficient
model (str, optional) – Regression model to fit to data.
Notes
model should be one of
ra fit ‘best’ robust model depending on whether or no spine width is less than slim
rs spine fit of Powell et al. (2020)
r2 robust model 2, a robust version of the Isoplot model 2
rx spine fit with expanded errors, for use if s is a little bit over slim (experimental feature only)
Data should be input as 1-dimensional ndarrays.
References
Powell, R., Green, E.C.R., Marillo Sialer, E., Woodhead, J., 2020. Robust isochron calculation. Geochronology 2, 325–342. https://doi.org/10.5194/gchron-2-325-2020
- pysoplot.regression.robust_model_2(x, y, xy_err=None, h=1.4, disp=True)
Fit a robust model 2 linear regression to 2-dimensional data.
Analytical errors discarded and takes the geometric mean slope of the spine y on x fit, and the reciprocal of the spine x on y fit, using data scatter in place of analytical weights. Similar to the Isoplot model 2 approach but with robust properties.
- Parameters
x (np.ndarray) – x values (as 1-dimensional array)
y (np.ndarray) – y values
xy_err (np.ndarray, optional) – 3 by n array with elements sx, sy, r_xy (i.e. analytical errors on x and y, and their correlation coefficient)
Notes
If xy_err is provided, nmad of the residuals will be checked against the median analytical uncertainty, and if larger, this value will be used to weight data points instead. Following Powell et al. (2002), this is to guard against a situation in which the scatter on the data is smaller than that expected with the analytical uncertainties.
Code by Roger Powell with style modifications by TP.
References
Powell, R., Hergt, J., Woodhead, J., 2002. Improving isochron calculations with robust statistics and the bootstrap. Chemical Geology 191–204.
- pysoplot.regression.siegel(x, y)
Median of pairwise median slopes algorithm of Siegel (1982).
Based on code by Roger Powell.
References
Siegel, A. F.: Robust regression Using repeated medians, Biometrika, 69, 242–244, 1982.
- pysoplot.regression.slim(n)
Upper bound of 95% confidence interval on s (spine width) for Spine linear regression algorithm. Derived from simulation of Gaussian distributed datasets. See Powell et al., (2020).
References
Powell, R., Green, E.C.R., Marillo Sialer, E., Woodhead, J., 2020. Robust isochron calculation. Geochronology 2, 325–342. https://doi.org/10.5194/gchron-2-325-2020
- pysoplot.regression.spine(x, sx, y, sy, rxy, h=1.4)
Iteratively re-weighted Huber line-fitting algorithm of Powell et al. (2020).
- Parameters
x (np.ndarray) – x values (as 1-dimensional array)
sx (np.ndarray) – analytical uncertainty on x at \(1\sigma\) abs.
y (np.ndarray) – y values
sy (np.ndaray) – analytical uncertainty on y at \(1\sigma\) abs.
rxy (np.ndarray) – x-y correlation coefficient
Notes
Updated to use normalised deltheta as convergence criteria for improved performance where slope and y-int are of very different order of magnitude.
Code by Roger Powell with style modifications by TP.
References
Powell, R., Green, E.C.R., Marillo Sialer, E., Woodhead, J., 2020. Robust isochron calculation. Geochronology 2, 325–342. https://doi.org/10.5194/gchron-2-325-2020
- pysoplot.regression.york(x, sx, y, sy, r_xy, sy_excess=None, model='1', itmax=50, rtol=1e-08, atol=1e-08, theta0=None)
Fit total weighted least squares regression using the algorithm of York (2004). Analytical errors on regression parameters are equivalent to the Maximum Likelihood approach of Titterington and Halliday (1979).
Can also fit a version of the Isoplot model 2 and model 3 (but see note) by adjusting th weightings. Model 2 code is based on McSaveney in Faure (1977). Model 3 code emulates Isoplot (Ludwig, 2012).
- Parameters
x (np.ndarray, 1-D) – x values
sx (np.ndarray, 1-D) – analytical uncertainty on x at \(1\sigma\) abs.
y (np.ndarray, 1-D) – y values
sy (np.ndaray, 1-D) – analytical uncertainty on y at \(1\sigma\) abs.
r_xy (np.ndarray, 1-D) – x-y correlation coefficient
model ({'1', '2', '3'}, optional) – model type to fit
Notes
Model 3 fits involve an iterative approach that typically requires calling this function several times. To perform a model 3 fit, call
pysoplot.regression.model_3()
instead.References
Faure, G., 1977. Appendix 1: Fitting of isochrons for dating by the Rb-Sr method, in: Principles of Isotope Geology. John Wiley and Sons, pp. 1–17.
Ludwig, K.R., 2012. Isoplot/Ex Version 3.75: A Geochronological Toolkit for Microsoft Excel, Special Publication 4. Berkeley Geochronology Center.
Titterington, D.M., Halliday, A.N., 1979. On the fitting of parallel isochrons and the method of maximum likelihood. Chemical Geology 26, 183–195.
York, D., Evensen, N.M., Martinez, M.L., De Basabe Delgado, J., 2004. Unified equations for the slope, intercept, and standard errors of the best straight line. American Journal of Physics 72, 367–375. https://doi.org/10.1119/1.1632486
transform
Transform data points, regression fits, etc.
- pysoplot.transform.dp_errors(dp, in_error_type, append_corr=True, row_wise=True, dim=2, tol=1e-10)
Convert data point uncertainties to 1 \(\sigma\) absolute.
- Parameters
in_error_type (str) – Error type and sigma level, use: {abs/per/rel}{1s/2s}
append_corr (bool, optional) – Assume error correlations are zero by appending a column of zeros to the end of n x 4 data point array.
row_wise (bool, optional) – True if each row is a data point and each column is a variable. False if vice versa.
dim (int) – 1 for univariate data (e.g. Pb*/U data), 2 for multivariate (2-D) data.
upb
Functions and routines for U-Pb geochronology.
Notes
See pysoplot.dqpb
for equivalent functions that account for
disequilibrium in the uranium-series decay chains.
- pysoplot.upb.concint_age(fit, method='L1980', diagram='tw', dc_errors=False)
Compute concordia-intercept age and uncertainty using various algorithms.
- Parameters
fit (
dict
) – Regression fit parameters.method ({'Powell', 'L1980', 'L2000'}, optional) – Algorithm used to compute age and age uncertainties.
- Returns
results – Age results.
- Return type
dict
Notes
The ‘Powell’ method focuses on lower intercept age only (< ~1 Ga), but works with all fit types. Age uncertainties are symmetric.
The ‘L1980’ method only works with classical fits. Age uncertainties are asymmetric.
The ‘L2000’ method works with all fit types and optionally accounts for decay constant errors, but requires 2 intercept age solutions.
- pysoplot.upb.isochron_age(fit, age_type='iso-206Pb', dc_errors=False, norm_isotope='204Pb')
Classical U-Pb isochron age and uncertainty.
- Parameters
fit (
dict
) – Regression fit parameters.age_type ({'iso-206Pb', 'iso-207Pb'}, optional) – Isochron type.
- Returns
results – Age results.
- Return type
dict
- pysoplot.upb.mc_concint(t, fit, trials=50000, diagram='tw', dc_errors=False, U_errors=False, intercept_plot=False, hist=False, xlim=(None, None), ylim=(None, None), env=False, age_ellipses=False, marker_max=None, marker_ages=(), auto_marker_ages=True, remove_overlaps=False, intercept_points=True, intercept_ellipse=False, negative_ages=True, age_prefix='Ma')
Compute Monte Carlo age uncertainties for equilibrium concordia intercept age.
- Parameters
t (float) – Calculated age (Ma).
fit (dict) – Linear regression fit parameters.
trials (int) – Number of Monte Carlo trials.
wtd_average
Weighted average algorithms for 1- and 2-dimensional data
- pysoplot.wtd_average.classical_wav(x, sx=None, V=None, method='ca')
Compute a classical 1-dimensional weighted average accounting for assigned uncertainties (and optionally uncertainty correlations).
Partly emulates the behaviour of Isoplot Ludwig (2012). If one-sided MSWD confidence limit is above lower threshold, uncertainty on weighted average is expanded according to data scatter (i.e., by sqrt(mswd)).
- Parameters
sx (np.ndarray, optional) – Analytical uncertainties at the 1 sigma level.
V (np.ndarray, optional) – Uncertainty covariance matrix.
method ({'ca'}) – Weighted average method to use (only one available currently)
Notes
Does not yet include external error component if MSWD confidence limit is above an upper threshold (as in the Isoplot MLE approach), although this feature may be added in future.
References
Ludwig, K.R., 2012. Isoplot/Ex Version 3.75: A Geochronological Toolkit for Microsoft Excel, Special Publication 4. Berkeley Geochronology Center.
- pysoplot.wtd_average.robust_wav(x, sx=None, V=None, method='ra')
Compute a robust 1-dimensional weighted average accounting for assigned uncertainties (and optionally uncertainty correlations).
- Parameters
x (np.ndarray) – Data points to average (as 1-d array).
sx (np.ndarray, optional) – Analytical uncertainties at the 1 sigma level (as 1-d array).
V (np.ndarray, optional) – Uncertainty covariance matrix.
method ({'ra'}) – Weighted average method to use (only one available currently)
Notes
Only implements spine robust weighted average at present.
- pysoplot.wtd_average.spine_wav_cor(x, V, xb0=None, maxiter=50, h=1.4, atol=1e-08, rtol=1e-08)
Compute a spine robust uncertainty weighted average accounting for uncertainty correlations amongst data points. Equivalent to a 1-dimensional version of the spine linear regression algorithm of Powell et al. (2020). Note, this reduces to classical statistics uncertainty weighted average for “well-behaved” datasets.
- Parameters
x (np.ndarray) – Data points to be averaged (1-dimensional array).
V (np.ndarray) – covariance matrix (2-dimensional array).
References
Powell, R., Green, E.C.R., Marillo Sialer, E., Woodhead, J., 2020. Robust isochron calculation. Geochronology 2, 325–342. https://doi.org/10.5194/gchron-2-325-2020
- pysoplot.wtd_average.wav(x, sx)
Compute a 1-dimensional uncertainty weighted mean without accounting for uncertainty covariances.
- pysoplot.wtd_average.wav_cor(x, V)
Compute a 1-dimensional uncertainty weighted mean, accounting for uncertainty correlations amongst individual data points.
Based on equations given in, e.g., Powell1 and Holland (1988], Lyons et al. (1988), and McLean et al. (2011).
- Parameters
x (np.ndarray) – Data points to be averaged (should be a 1-dimensional array).
V (np.ndarray) – covariance matrix (should be a 2-dimensional array).
Notes
See cautionary note in Lyons et al. (1988):
“One point to beware of in practice is the situation in which the individual weights become numerically very large.” … this happens when [errors are very close to each other and highly correlated]. If we have slightly mis-estimated the elements of the error matrix, or if our measurements are slightly biassed, the effect of the large weights with different signs can be to drive our solution far away from the correct value.”
References
Lyons, L., Gibaut, D., Clifford, P., 1988. How to combine correlated estimates of a single physical quantity. Nuclear Instruments and Methods in Physics Research Section A: Accelerators, Spectrometers, Detectors and Associated Equipment 270, 110–117.
McLean, N M, J F Bowring, and S A Bowring, 2011, An Algorithm for U-Pb Isotope Dilution Data Reduction and Uncertainty Propagation Geochemistry, Geophysics, Geosystems 12, no. 6. https://doi.org/10.1029/2010GC003478.
Powell, R., Holland, T.J.B., 1988. An internally consistent dataset with uncertainties and correlations: 3. Applications to geobarometry, worked examples and a computer program. J Metamorph Geol 6, 173–204. https://doi.org/10.1111/j.1525-1314.1988.tb00415.x
ludwig
Disequilibrium U-Pb equations based on Ludwig (1977).
- pysoplot.ludwig.bateman(Lam, series='238U')
Return Bateman coefficients for the 238U or 235U decay series as defined in Ludwig (1977).
- Lamarray-like
array of relevant decay constants (see below)
Notes
References
Ludwig, K.R., 1977. Effect of initial radioactive-daughter disequilibrium on U-Pb isotope apparent ages of young minerals. Journal of Research of the US Geological Survey 5, 663–667.
- pysoplot.ludwig.f(t, A, Lam, coef, init=(True, True))
206Pb*/238U ratio (where * denotes radiogenic Pb) as a function of t and activity ratio values following Ludwig (1977). Note there is a small typo in the original Ludwig article.
- Parameters
t (float or array-like) – Age (Ma)
A (array-like) – One-dimensional array of activity ratio values arranged as follows [234U/238U], [230Th/238U], [226Ra/238U],
Lam (array-like) – One-dimensional array of decay constants (Ma^-1) arragend as follows [lam238, lam234, lam230, lam226]
coef (array-like) – One-dimensional array of Bateman coefficients [c1, c2, c3, c4, c5, p1, p2, p3, h1, h2]
init (array-like, optional) – Two element array of bools. First element is True if [234U/238U] is an initial value. Second element is True if [230Th/238U] is an initial value.
References
Ludwig, K.R., 1977. Effect of initial radioactive-daughter disequilibrium on U-Pb isotope apparent ages of young minerals. Journal of Research of the US Geological Survey 5, 663–667.
- pysoplot.ludwig.g(t, a231_235_i, Lam, coef)
207Pb*/235U ratio (where * denotes radiogenic Pb) as a function of t and activity ratio values following Ludwig (1977).
- Parameters
t (float or array-like) – Age (Ma)
a231_235_i (float or array-like) – [231Pa/235U] activity ratio.
Lam (array-like) – One-dimensional array of decay constants (Ma^-1) arragend as follows [lam235, lam231]
coef (array-like) – One-dimensional array of Bateman coefficients [d1, d2]
References
Ludwig, K.R., 1977. Effect of initial radioactive-daughter disequilibrium on U-Pb isotope apparent ages of young minerals. Journal of Research of the US Geological Survey 5, 663–667.