pyrolite.util.lambdas
- pyrolite.util.lambdas.calc_lambdas(df, params=None, degree=4, exclude=[], algorithm='ONeill', anomalies=[], fit_tetrads=False, sigmas=None, add_uncertainties=False, add_X2=False, **kwargs)[source]
Parameterises values based on linear combination of orthogonal polynomials over a given set of values for independent variable x 1 . This function expects to recieve data already normalised and log-transformed.
- Parameters
df (
pd.DataFrame
) – Dataframe containing REE Data.params (
list
|str
) – Pre-computed parameters for the orthogonal polynomials (a list of tuples). Optionally specified, otherwise defaults the parameterisation as in O’Neill (2016). 1 If a string is supplied,"O'Neill (2016)"
or similar will give the original defaults, while"full"
will use all of the REE (including Eu) as a basis for the orthogonal polynomials.degree (
int
) – Degree of orthogonal polynomial fit.exclude (
list
) – REE to exclude from the fit.algorithm (
str
) – Algorithm to use for fitting the orthogonal polynomials.anomalies (
list
) – List of relative anomalies to append to the dataframe.fit_tetrads (
bool
) – Whether to fit tetrad functions in addition to orthogonal polynomial functions. This will force the use of the optimization algorithm.sigmas (
float
|numpy.ndarray
) – Single value or 1D array of observed value uncertainties.add_uncertainties (
bool
) – Whether to append estimated parameter uncertainties to the dataframe.add_X2 (
bool
) – Whether to append the chi-squared values (χ2) to the dataframe.- Return type
pd.DataFrame
References
- 1(1,2)
O’Neill HSC (2016) The Smoothness and Shapes of Chondrite-normalized Rare Earth Element Patterns in Basalts. J Petrology 57:1463–1508. doi: 10.1093/petrology/egw047
pyrolite.util.lambdas.params
Functions to generate parameters for the construction of orthogonal polynomials which are used to fit REE patterns.
- pyrolite.util.lambdas.params.orthogonal_polynomial_constants(xs, degree=3, rounding=None, tol=1e-14)[source]
Finds the parameters \((\beta_0), (\gamma_0, \gamma_1), (\delta_0, \delta_1, \delta_2)\) etc. for constructing orthogonal polynomial functions f(x) over a fixed set of values of independent variable x. Used for obtaining lambda values for dimensional reduction of REE data 2.
- Parameters
xs (
numpy.ndarray
) – Indexes over which to generate the orthogonal polynomials.degree (
int
) – Maximum polynomial degree. E.g. 2 will generate constant, linear, and quadratic polynomial components.tol (
float
) – Convergence tolerance for solver.rounding (
int
) – Precision for the orthogonal polynomial coefficents.- Returns
List of tuples corresponding to coefficients for each of the polynomial components. I.e the first tuple will be empty, the second will contain a single coefficient etc.
- Return type
Notes
Parameters are used to construct orthogonal polymomials of the general form:
\[\begin{split}f(x) &= a_0 \\ &+ a_1 * (x - \beta) \\ &+ a_2 * (x - \gamma_0) * (x - \gamma_1) \\ &+ a_3 * (x - \delta_0) * (x - \delta_1) * (x - \delta_2) \\\end{split}\]See also
References
- 2
O’Neill HSC (2016) The Smoothness and Shapes of Chondrite-normalized Rare Earth Element Patterns in Basalts. J Petrology 57:1463–1508. doi: 10.1093/petrology/egw047
- pyrolite.util.lambdas.params.parse_sigmas(size, sigmas=None)[source]
Disambigaute a value or set of sigmas for a dataset for use in lambda-fitting algorithms.
- Parameters
sigmas (
float
|numpy.ndarray
) – 2D array of REE uncertainties. Values as fractional uncertaintes (i.e. \(\sigma_{REE} / REE\)).- Returns
sigmas – 1D array of sigmas (\(\sigma_{REE} / REE\)).
- Return type
Notes
Note that the y-array is passed here only to infer the shape which should be assumed by the uncertainty array. Through propagation of uncertainties, the uncertainty on the natural logarithm of the normalised REE values are equivalent to \(\sigma_{REE} / REE\) where the uncertainty in the reference composition is assumed to be zero. Thus, standard deviations of 1% in REE will result in \(\sigma=0.01\) for the log-transformed REE. If no sigmas are provided, 1% uncertainty will be assumed and an array of 0.01 will be returned.
pyrolite.util.lambdas.eval
Generation and evalutation of orthogonal polynomial and tetrad functions from sets of parameters (the sequence of polymomial roots and tetrad centres and widths).
- pyrolite.util.lambdas.eval.lambda_poly(x, ps)[source]
Evaluate polynomial lambda_n(x) given a tuple of parameters ps with length equal to the polynomial degree.
- Parameters
x (
numpy.ndarray
) – X values to calculate the function at.ps (
tuple
) – Parameter set tuple. E.g. parameters (a, b) from \(f(x) = (x-a)(x-b)\).- Return type
- pyrolite.util.lambdas.eval.tetrad(x, centre, width)[source]
Evaluate \(f(z)\) describing a tetrad with specified centre and width.
- pyrolite.util.lambdas.eval.get_lambda_poly_function(lambdas: ndarray, params=None, radii=None, degree=5)[source]
Expansion of lambda parameters back to the original space. Returns a function which evaluates the sum of the orthogonal polynomials at given x values.
- Parameters
lambdas (
numpy.ndarray
) – Lambda values to weight combination of polynomials.params (
list
(tuple
)) – Parameters for the orthogonal polynomial decomposition.radii (
numpy.ndarray
) – Radii values used to construct the lambda values. 3degree (
int
) – Degree of the orthogonal polynomial decomposition. 3Notes
pyrolite.util.lambdas.oneill
Linear algebra methods for fitting a series of orthogonal polynomial functions to REE patterns.
- pyrolite.util.lambdas.oneill.get_polynomial_matrix(radii, params=None)[source]
Create the matrix A with polynomial components across the columns, and increasing order down the rows.
- Parameters
radii (
list
,numpy.ndarray
) – Radii at which to evaluate the orthogonal polynomial.params (
tuple
) – Tuple of constants for the orthogonal polynomial.- Return type
See also
- pyrolite.util.lambdas.oneill.lambdas_ONeill2016(df, radii, params=None, sigmas=None, add_X2=False, add_uncertainties=False, **kwargs)[source]
Implementation of the original algorithm. 4
- Parameters
df (
pandas.DataFrame
|pandas.Series
) – Dataframe of REE data, with sample analyses organised by row.radii (
list
,numpy.ndarray
) – Radii at which to evaluate the orthogonal polynomial.params (
tuple
) – Tuple of constants for the orthogonal polynomial.sigmas (
float
|numpy.ndarray
) – Single value or 1D array of normalised observed value uncertainties (\(\sigma_{REE} / REE\)).add_X2 (
bool
) – Whether to append the chi-squared values (χ2) to the dataframe/series.add_uncertainties (
bool
) – Append parameter standard errors to the dataframe/series.- Return type
References
- 4
O’Neill HSC (2016) The Smoothness and Shapes of Chondrite-normalized Rare Earth Element Patterns in Basalts. J Petrology 57:1463–1508. doi: 10.1093/petrology/egw047
pyrolite.util.lambdas.opt
Functions for optimization-based REE profile fitting and parameter uncertainty estimation.
- pyrolite.util.lambdas.opt.pcov_from_jac(jac)[source]
Extract a covariance matrix from a Jacobian matrix returned from
scipy.optimize
functions.
- Parameters
jac (
numpy.ndarray
) – Jacobian array.- Returns
pcov – Square covariance array; this hasn’t yet been scaled by residuals.
- Return type
- pyrolite.util.lambdas.opt.linear_fit_components(y, x0, func_components, sigmas=None)[source]
Fit a weighted sum of function components using linear algebra.
- Parameters
y (
numpy.ndarray
) – Array of target values to fit.x0 (
numpy.ndarray
) – Starting guess for the function weights.func_components (
list
(numpy.ndarray
)) – List of arrays representing static/evaluated function components.sigmas (
float
|numpy.ndarray
) – Single value or 1D array of normalised observed value uncertainties (\(\sigma_{REE} / REE\)).- Returns
B, s, χ2 – Arrays for the optimized parameter values (B; (n, d)), parameter uncertaintes (s, 1σ; (n, d)) and chi-chi_squared (χ2; (n, 1)).
- Return type
- pyrolite.util.lambdas.opt.optimize_fit_components(y, x0, func_components, residuals_function=<function _residuals_func>, sigmas=None)[source]
Fit a weighted sum of function components using
scipy.optimize.least_squares()
.
- Parameters
y (
numpy.ndarray
) – Array of target values to fit.x0 (
numpy.ndarray
) – Starting guess for the function weights.func_components (
list
(numpy.ndarray
)) – List of arrays representing static/evaluated function components.redsiduals_function (callable) – Callable funciton to compute residuals which accepts ordered arguments for weights, target values and function components.
sigmas (
float
|numpy.ndarray
) – Single value or 1D array of normalised observed value uncertainties (\(\sigma_{REE} / REE\)).- Returns
B, s, χ2 – Arrays for the optimized parameter values (B; (n, d)), parameter uncertaintes (s, 1σ; (n, d)) and chi-chi_squared (χ2; (n, 1)).
- Return type
- pyrolite.util.lambdas.opt.lambdas_optimize(df, radii, params=None, fit_tetrads=False, tetrad_params=None, fit_method='opt', sigmas=None, add_uncertainties=False, add_X2=False, **kwargs)[source]
Parameterises values based on linear combination of orthogonal polynomials over a given set of values for independent variable x. 5
- Parameters
df (
pandas.DataFrame
| :class:`pandas.Series) – Target data to fit. For geochemical data, this is typically normalised so we can fit a smooth function.radii (
list
,numpy.ndarray
) – Radii at which to evaluate the orthogonal polynomial.params (
list
,None
) – Orthogonal polynomial coefficients (seeorthogonal_polynomial_constants()
).fit_tetrads (
bool
) – Whether to also fit the patterns for tetrads.tetrad_params (
list
) – List of parameter sets for tetrad functions.fit_method (
str
) – Which fit method to use:"optimization"
or"linear"
.sigmas (
float
|numpy.ndarray
) – Single value or 1D array of observed value uncertainties.add_uncertainties (
bool
) – Whether to append estimated parameter uncertainties to the dataframe/series.add_X2 (
bool
) – Whether to append the chi-squared values (χ2) to the dataframe/series.- Returns
Optimial results for weights of orthogonal polymomial regression (lambdas).
- Return type
References
- 5
O’Neill HSC (2016) The Smoothness and Shapes of Chondrite-normalized Rare Earth Element Patterns in Basalts. J Petrology 57:1463–1508. doi: 10.1093/petrology/egw047
pyrolite.util.lambdas.plot
Functions for the visualisation of reconstructed and deconstructed parameterised REE profiles based on parameterisations using ‘lambdas’ (and tetrad-equivalent weights ‘taus’).
- pyrolite.util.lambdas.plot.plot_lambdas_components(lambdas, params=None, ax=None, **kwargs)[source]
Plot a decomposed orthogonal polynomial from a single set of lambda coefficients.
- Parameters
lambdas – 1D array of lambdas.
params (
list
) – List of orthongonal polynomial parameters, if defaults are not used.ax (
matplotlib.axes.Axes
) – Optionally specified axes to plot on.index (
str
) – Index to use for the plot (one of"index", "radii", "z"
).- Return type
- pyrolite.util.lambdas.plot.plot_tetrads_components(taus, tetrad_params=None, ax=None, index='radii', logy=True, drop0=True, **kwargs)[source]
Individually plot the four tetrad components for one set of $ au$s.
- Parameters
taus (
numpy.ndarray
) – 1D array of $ au$ tetrad function coefficients.tetrad_params (
list
) – List of tetrad parameters, if defaults are not used.ax (
matplotlib.axes.Axes
) – Optionally specified axes to plot on.index (
str
) – Index to use for the plot (one of"index", "radii", "z"
).logy (
bool
) – Whether to log-scale the y-axis.drop0 (
bool
) – Whether to remove zeroes from the outputs such that individual tetrad functions are shown only within their respective bounds (and not across the entire REE, where their effective values are zero).
- pyrolite.util.lambdas.plot.plot_profiles(coefficients, tetrads=False, params=None, tetrad_params=None, ax=None, index='radii', logy=False, **kwargs)[source]
Plot the reconstructed REE profiles of a 2D dataset of coefficients ($lambda$s, and optionally $tau$s).
- Parameters
coefficients (
numpy.ndarray
) – 2D array of $lambda$ orthogonal polynomial coefficients, and optionally including $tau$ tetrad function coefficients in the last four columns (wheretetrads=True
).tetrads (
bool
) – Whether the coefficient array contains tetrad coefficients ($tau$s).params (
list
) – List of orthongonal polynomial parameters, if defaults are not used.tetrad_params (
list
) – List of tetrad parameters, if defaults are not used.ax (
matplotlib.axes.Axes
) – Optionally specified axes to plot on.index (
str
) – Index to use for the plot (one of"index", "radii", "z"
).logy (
bool
) – Whether to log-scale the y-axis.- Return type
pyrolite.util.lambdas.transform
Functions for transforming ionic radii to and from atomic number for the visualisation of REE patterns.
- pyrolite.util.lambdas.transform.REE_z_to_radii(z, fit=None, degree=7, **kwargs)[source]
Estimate the ionic radii which would be approximated by a given atomic number based on a provided (or calcuated) fit for the Rare Earth Elements.
- Parameters
z (
float
|list
|numpy.ndarray
) – Atomic nubmers to be converted.fit (callable) – Callable function optionally specified; if not specified it will be calculated from Shannon Radii.
degree (
int
) – Degree of the polynomial fit between atomic number and radii.- Returns
r – Approximate atomic nubmers for given radii.
- Return type
- pyrolite.util.lambdas.transform.REE_radii_to_z(r, fit=None, degree=7, **kwargs)[source]
Estimate the atomic number which would be approximated by a given ionic radii based on a provided (or calcuated) fit for the Rare Earth Elements.
- Parameters
r (
float
|list
|numpy.ndarray
) – Radii to be converted.fit (callable) – Callable function optionally specified; if not specified it will be calculated from Shannon Radii.
degree (
int
) – Degree of the polynomial fit between radii and atomic number.- Returns
z – Approximate atomic numbers for given radii.
- Return type