Title: | Statistical Toolbox for Sedimentary Provenance Analysis |
---|---|
Description: | Bundles a number of established statistical methods to facilitate the visual interpretation of large datasets in sedimentary geology. Includes functionality for adaptive kernel density estimation, principal component analysis, correspondence analysis, multidimensional scaling, generalised procrustes analysis and individual differences scaling using a variety of dissimilarity measures. Univariate provenance proxies, such as single-grain ages or (isotopic) compositions are compared with the Kolmogorov-Smirnov, Kuiper, Wasserstein-2 or Sircombe-Hazelton L2 distances. Categorical provenance proxies such as chemical compositions are compared with the Aitchison and Bray-Curtis distances,and count data with the chi-square distance. Varietal data can either be converted to one or more distributional datasets, or directly compared using the multivariate Wasserstein distance. Also included are tools to plot compositional and count data on ternary diagrams and point-counting data on radial plots, to calculate the sample size required for specified levels of statistical precision, and to assess the effects of hydraulic sorting on detrital compositions. Includes an intuitive query-based user interface for users who are not proficient in R. |
Authors: | Pieter Vermeesch [aut, cre] |
Maintainer: | Pieter Vermeesch <[email protected]> |
License: | GPL-2 |
Version: | 4.4 |
Built: | 2025-02-10 06:29:59 UTC |
Source: | https://github.com/pvermees/provenance |
Calculates Aitchison's additive logratio transformation for a
dataset of class compositional
or a compositional data
matrix.
ALR(x, ...) ## Default S3 method: ALR(x, inverse = FALSE, ...) ## S3 method for class 'compositional' ALR(x, ...)
ALR(x, ...) ## Default S3 method: ALR(x, inverse = FALSE, ...) ## S3 method for class 'compositional' ALR(x, ...)
x |
an object of class |
... |
optional arguments |
inverse |
perform the inverse inverse logratio transformation? |
a matrix of ALR coordinates OR an object of class
compositional
(if inverse=TRUE
).
# logratio plot of trace element concentrations: data(Namib) alr <- ALR(Namib$Trace) pairs(alr[,1:5]) title('log(X/Pb)')
# logratio plot of trace element concentrations: data(Namib) alr <- ALR(Namib$Trace) pairs(alr[,1:5]) title('log(X/Pb)')
Adds several components of a composition together into a single component
amalgamate(x, ...) ## Default S3 method: amalgamate(x, ...) ## S3 method for class 'compositional' amalgamate(x, ...) ## S3 method for class 'counts' amalgamate(x, ...) ## S3 method for class 'SRDcorrected' amalgamate(x, ...) ## S3 method for class 'varietal' amalgamate(x, ...)
amalgamate(x, ...) ## Default S3 method: amalgamate(x, ...) ## S3 method for class 'compositional' amalgamate(x, ...) ## S3 method for class 'counts' amalgamate(x, ...) ## S3 method for class 'SRDcorrected' amalgamate(x, ...) ## S3 method for class 'varietal' amalgamate(x, ...)
x |
a compositional dataset |
... |
a series of new labels assigned to strings or vectors of strings denoting the components that need amalgamating |
an object of the same class as X
with fewer
components
data(Namib) HMcomponents <- c("zr","tm","rt","TiOx","sph","ap","ep", "gt","st","amp","cpx","opx") am <- amalgamate(Namib$PTHM,feldspars=c("KF","P"), lithics=c("Lm","Lv","Ls"),heavies=HMcomponents) plot(ternary(am))
data(Namib) HMcomponents <- c("zr","tm","rt","TiOx","sph","ap","ep", "gt","st","amp","cpx","opx") am <- amalgamate(Namib$PTHM,feldspars=c("KF","P"), lithics=c("Lm","Lv","Ls"),heavies=HMcomponents) plot(ternary(am))
acomp
objectConvert an object of class compositional
to an object of
class acomp
for use in the compositions
package
as.acomp(x)
as.acomp(x)
x |
an object of class |
a data.frame
data(Namib) qfl <- ternary(Namib$PT,c('Q'),c('KF','P'),c('Lm','Lv','Ls')) plot(qfl,type="QFL.dickinson") qfl.acomp <- as.acomp(qfl) ## uncomment the next two lines to plot an error ## ellipse using the 'compositions' package: # library(compositions) # ellipses(mean(qfl.acomp),var(qfl.acomp),r=2)
data(Namib) qfl <- ternary(Namib$PT,c('Q'),c('KF','P'),c('Lm','Lv','Ls')) plot(qfl,type="QFL.dickinson") qfl.acomp <- as.acomp(qfl) ## uncomment the next two lines to plot an error ## ellipse using the 'compositions' package: # library(compositions) # ellipses(mean(qfl.acomp),var(qfl.acomp),r=2)
compositional
objectConvert an object of class matrix
, data.frame
or
acomp
to an object of class compositional
as.compositional(x, method = NULL, colmap = "rainbow")
as.compositional(x, method = NULL, colmap = "rainbow")
x |
an object of class |
method |
dissimilarity measure, either |
colmap |
the colour map to be used in pie charts. |
an object of class compositional
data(Namib) PT.acomp <- as.acomp(Namib$PT) PT.compositional <- as.compositional(PT.acomp) print(Namib$PT$x - PT.compositional$x) ## uncomment the following lines for an illustration of using this ## function to integrate 'provenance' with 'compositions' # library(compositions) # data(Glacial) # a.glac <- acomp(Glacial) # c.glac <- as.compositional(a.glac) # summaryplot(c.glac,ncol=8)
data(Namib) PT.acomp <- as.acomp(Namib$PT) PT.compositional <- as.compositional(PT.acomp) print(Namib$PT$x - PT.compositional$x) ## uncomment the following lines for an illustration of using this ## function to integrate 'provenance' with 'compositions' # library(compositions) # data(Glacial) # a.glac <- acomp(Glacial) # c.glac <- as.compositional(a.glac) # summaryplot(c.glac,ncol=8)
counts
objectConvert an object of class matrix
or data.frame
to an
object of class counts
as.counts(x, method = "chisq", colmap = "rainbow")
as.counts(x, method = "chisq", colmap = "rainbow")
x |
an object of class |
method |
either |
colmap |
the colour map to be used in pie charts. |
an object of class counts
X <- matrix(c(0,100,0,30,11,2,94,36,0),nrow=3,ncol=3) rownames(X) <- 1:3 colnames(X) <- c('a','b','c') comp <- as.counts(X) d <- diss(comp)
X <- matrix(c(0,100,0,30,11,2,94,36,0),nrow=3,ncol=3) rownames(X) <- 1:3 colnames(X) <- c('a','b','c') comp <- as.counts(X) d <- diss(comp)
data.frame
objectConvert an object of class compositional
to a
data.frame
for use in the robCompositions
package
## S3 method for class 'compositional' as.data.frame(x, ...) ## S3 method for class 'counts' as.data.frame(x, ...)
## S3 method for class 'compositional' as.data.frame(x, ...) ## S3 method for class 'counts' as.data.frame(x, ...)
x |
an object of class |
... |
optional arguments to be passed on to the generic function |
a data.frame
data(Namib) Major.frame <- as.data.frame(Namib$Major) ## uncomment the next two lines to plot an error ## ellipse using the robCompositions package: # library(robCompositions) # plot(pcaCoDa(Major.frame))
data(Namib) Major.frame <- as.data.frame(Namib$Major) ## uncomment the next two lines to plot an error ## ellipse using the robCompositions package: # library(robCompositions) # plot(pcaCoDa(Major.frame))
varietal
objectConvert an object of class matrix
or data.frame
to an
object of class varietal
as.varietal(x, snames = NULL, method = "KS")
as.varietal(x, snames = NULL, method = "KS")
x |
an object of class |
snames |
either a vector of sample names, an integer marking
the length of the sample name prefix, or |
method |
either |
an object of class varietal
fn <- system.file("SNSM/Ttn_chem.csv",package="provenance") ap1 <- read.csv(fn) ap2 <- as.varietal(x=ap1,snames=3)
fn <- system.file("SNSM/Ttn_chem.csv",package="provenance") ap1 <- read.csv(fn) ap2 <- as.varietal(x=ap1,snames=3)
Uses the diffusion algorithm of Botev (2011) to calculate the bandwidth for kernel density estimation
botev(x)
botev(x)
x |
a vector of ordinal data |
a scalar value with the optimal bandwidth
Zdravko Botev
Botev, Z. I., J. F. Grotowski, and D. P. Kroese. "Kernel density estimation via diffusion." The Annals of Statistics 38.5 (2010): 2916-2957.
fname <- system.file("Namib/DZ.csv",package="provenance") bw <- botev(read.distributional(fname)$x$N1) print(bw)
fname <- system.file("Namib/DZ.csv",package="provenance") bw <- botev(read.distributional(fname)$x$N1) print(bw)
Calculates the Bray-Curtis dissimilarity between two samples
bray.diss(x, ...) ## Default S3 method: bray.diss(x, y, ...) ## S3 method for class 'compositional' bray.diss(x, ...)
bray.diss(x, ...) ## Default S3 method: bray.diss(x, y, ...) ## S3 method for class 'compositional' bray.diss(x, ...)
x |
a vector containing the first compositional sample |
... |
optional arguments |
y |
a vector of |
a scalar value
data(Namib) print(bray.diss(Namib$HM$x["N1",],Namib$HM$x["N2",]))
data(Namib) print(bray.diss(Namib$HM$x["N1",],Namib$HM$x["N2",]))
Performs Correspondence Analysis of point-counting data
CA(x, nf = 2, ...)
CA(x, nf = 2, ...)
x |
an object of class |
nf |
number of correspondence factors (dimensions) |
... |
optional arguments to the |
an object of classes CA
, which is synonymous to the
MASS package's correspondence
class.
data(Namib) plot(CA(Namib$PT))
data(Namib) plot(CA(Namib$PT))
Computes the logratio mean composition of a continuous mixture of point-counting data.
## S3 method for class 'counts' central(x, ...)
## S3 method for class 'counts' central(x, ...)
x |
an object of class |
... |
optional arguments |
The central composition assumes that the observed point-counting distribution is the combination of two sources of scatter: counting uncertainty and true geological dispersion.
an [5 x n]
matrix with n
being the number
of categories and the rows containing:
the ‘central’ composition.
the standard error for the central composition.
the overdispersion parameter, i.e. the coefficient of
variation of the underlying logistic normal
distribution. central
computes a continuous
mixture model for each component (column)
separately. Covariance terms are not reported.
the lower limit of a ‘1 sigma’ region for theta
.
the upper limit of a ‘1 sigma’ region for theta
.
the mean square of the weighted deviates, a.k.a. reduced chi-square statistic.
the p-value for age homogeneity
Calculates Aitchison's centered logratio transformation for a
dataset of class compositional
or a compositional data
matrix.
CLR(x, ...) ## Default S3 method: CLR(x, inverse = FALSE, ...) ## S3 method for class 'compositional' CLR(x, ...)
CLR(x, ...) ## Default S3 method: CLR(x, inverse = FALSE, ...) ## S3 method for class 'compositional' CLR(x, ...)
x |
an object of class |
... |
optional arguments |
inverse |
perform the inverse inverse logratio transformation? |
a matrix of CLR coordinates OR an object of class
compositional
(if inverse=TRUE
)
# The following code shows that applying provenance's PCA function # to compositional data is equivalent to applying R's built-in # princomp function to the CLR transformed data. data(Namib) plot(PCA(Namib$Major)) dev.new() clrdat <- CLR(Namib$Major) biplot(princomp(clrdat))
# The following code shows that applying provenance's PCA function # to compositional data is equivalent to applying R's built-in # princomp function to the CLR transformed data. data(Namib) plot(PCA(Namib$Major)) dev.new() clrdat <- CLR(Namib$Major) biplot(princomp(clrdat))
Lumps all single grain analyses of several samples together under a new name
combine(x, ...)
combine(x, ...)
x |
a distributional dataset |
... |
a series of new labels assigned to strings or vectors of strings denoting the samples that need amalgamating |
a distributional data object with fewer samples than
x
data(Namib) combined <- combine(Namib$DZ, east=c('N3','N4','N5','N6','N7','N8','N9','N10'), west=c('N1','N2','N11','N12','T8','T13')) summaryplot(KDEs(combined))
data(Namib) combined <- combine(Namib$DZ, east=c('N3','N4','N5','N6','N7','N8','N9','N10'), west=c('N1','N2','N11','N12','T8','T13')) summaryplot(KDEs(combined))
List of rock and mineral densities using the following abbreviations: Q (quartz), KF (K-feldspar), P (plagioclase), F (feldspar), Lvf (felsic/porfiritic volcanic rock fragments), Lvm (microlithic / porfiritic / trachitic volcanic rock fragments), Lcc (calcite), Lcd (dolomite), Lp (marl), Lch (chert), Lms (argillaceous / micaceous rock fragments), Lmv (metavolcanics), Lmf (metasediments), Lmb (metabasites), Lv (volcanic rock fragments), Lc (carbonates), Ls (sedimentary rock fragments), Lm (metamorphic rock fragments), Lu (serpentinite), mica, opaques, FeOx (Fe-oxides), turbids, zr (zircon), tm (tourmaline), rt (rutile), TiOx (Ti-oxides), sph (titanite), ap (apatite), mon (monazite), oth (other minerals), ep (epidote), othLgM (prehnite + pumpellyite + lawsonite + carpholite), gt (garnet), ctd (chloritoid), st (staurolite), and (andalusite), ky (kyanite), sil (sillimanite), amp (amphibole), px (pyroxene), cpx (clinopyroxene), opx (orthopyroxene), ol (olivine), spinel and othHM (other heavy minerals).
Alberto Resentini and Pieter Vermeesch
Resentini, A, Malusa M G and Garzanti, E. "MinSORTING: An Excel worksheet for modelling mineral grain-size distribution in sediments, with application to detrital geochronology and provenance studies." Computers & Geosciences 59 (2013): 90-97.
Garzanti, E, Ando, S and Vezzoli, G. "Settling equivalence of detrital minerals and grain-size dependence of sediment composition." Earth and Planetary Science Letters 273.1 (2008): 138-151.
restore, minsorting
N8 <- subset(Namib$HM,select="N8") distribution <- minsorting(N8,densities,phi=2,sigmaphi=1,medium="air",by=0.05) plot(distribution)
N8 <- subset(Namib$HM,select="N8") distribution <- minsorting(N8,densities,phi=2,sigmaphi=1,medium="air",by=0.05) plot(distribution)
distributional
, compositional
, counts
or
varietal
Calculate the dissimilarity matrix between two datasets of class
distributional
or compositional
using the Kolmogorov-Smirnov,
Sircombe-Hazelton, Aitchison or Bray-Curtis distance
## S3 method for class 'distributional' diss(x, method = NULL, log = FALSE, verbose = FALSE, ...) ## S3 method for class 'compositional' diss(x, method = NULL, ...) ## S3 method for class 'counts' diss(x, method = NULL, ...) ## S3 method for class 'varietal' diss(x, method = NULL, ...)
## S3 method for class 'distributional' diss(x, method = NULL, log = FALSE, verbose = FALSE, ...) ## S3 method for class 'compositional' diss(x, method = NULL, ...) ## S3 method for class 'counts' diss(x, method = NULL, ...) ## S3 method for class 'varietal' diss(x, method = NULL, ...)
x |
an object of class |
method |
if if if if |
log |
logical. If |
verbose |
logical. If |
... |
optional arguments |
"KS"
stands for the Kolmogorov-Smirnov statistic,
"W2_1D"
for the 1-dimensional Wasserstein-2 distance,
"Kuiper"
for the Kuiper statistic, "SH"
for the
Sircombe-Hazelton distance, "aitchison"
for the
Aitchison logratio distance, "bray"
for the Bray-Curtis
distance, "chisq"
for the Chi-square distance, and "W2"
for the 2-dimensional Wasserstein-2 distance.
an object of class diss
KS.diss bray.diss SH.diss Wasserstein.diss Kuiper.diss
data(Namib) print(round(100*diss(Namib$DZ)))
data(Namib) print(round(100*diss(Namib$DZ)))
A compositional dataset comprising the mineralogical compositions
of the following end-members: undissected_magmatic_arc
,
dissected_magmatic_arc
, ophiolite
,
recycled_clastic
,
undissected_continental_block
,
transitional_continental_block
,
dissected_continental_block
,
subcreted_axial_belt
and
subducted_axial_belt
Alberto Resentini and Pieter Vermeesch
Resentini, A, Malusa M G and Garzanti, E. "MinSORTING: An Excel worksheet for modelling mineral grain-size distribution in sediments, with application to detrital geochronology and provenance studies." Computers & Geosciences 59 (2013): 90-97.
Garzanti, E, Ando, S and Vezzoli, G. "Settling equivalence of detrital minerals and grain-size dependence of sediment composition." Earth and Planetary Science Letters 273.1 (2008): 138-151.
minsorting
ophiolite <- subset(endmembers,select="ophiolite") plot(minsorting(ophiolite,densities,by=0.05))
ophiolite <- subset(endmembers,select="ophiolite") plot(minsorting(ophiolite,densities,by=0.05))
For a given sample size, returns the largest fraction which has been sampled with (1-p) x 100 % likelihood.
get.f(n, p = 0.05)
get.f(n, p = 0.05)
n |
the number of grains in the detrital sample |
p |
the required level of confidence |
the largest fraction that is sampled with at least (1-p) x 100% certainty
Vermeesch, Pieter. "How many grains are needed for a provenance study?" Earth and Planetary Science Letters 224.3 (2004): 441-451.
print(get.f(60)) print(get.f(117))
print(get.f(60)) print(get.f(117))
Returns the number of grains that need to be analysed to decrease the likelihood of missing any fraction greater than a given size below a given level.
get.n(p = 0.05, f = 0.05)
get.n(p = 0.05, f = 0.05)
p |
the probability that all n grains in the sample have missed
at least one fraction of size |
f |
the size of the smallest resolvable fraction (0<f<1) |
n |
the number of grains in the sample |
the number of grains needed to reduce the chance of missing
at least one fraction f of the total population to less than p
Vermeesch, Pieter. "How many grains are needed for a provenance study?." Earth and Planetary Science Letters 224.3 (2004): 441-451.
# number of grains required to be 99% that no fraction greater than 5% was missed: print(get.n(0.01)) # number of grains required to be 90% that no fraction greater than 10% was missed: print(get.n(p=0.1,f=0.1))
# number of grains required to be 99% that no fraction greater than 5% was missed: print(get.n(0.01)) # number of grains required to be 90% that no fraction greater than 10% was missed: print(get.n(p=0.1,f=0.1))
For a given sample size, returns the likelihood of missing any fraction greater than a given size
get.p(n, f = 0.05)
get.p(n, f = 0.05)
n |
the number of grains in the detrital sample |
f |
the size of the smallest resolvable fraction
(0< |
the probability that all n
grains in the sample have
missed at least one fraction of size f
Vermeesch, Pieter. "How many grains are needed for a provenance study?." Earth and Planetary Science Letters 224.3 (2004): 441-451.
print(get.p(60)) print(get.p(117))
print(get.p(60)) print(get.p(117))
Given a number of (2D) configurations, this function uses a combination of transformations (reflections, rotations, translations and scaling) to find a ‘consensus’ configuration which best matches all the component configurations in a least-squares sense.
GPA(X, scale = TRUE)
GPA(X, scale = TRUE)
X |
a list of dissimilarity matrices |
scale |
boolean flag indicating if the transformation should include the scaling operation |
a two column vector with the coordinates of the group configuration
procrustes
Performs 3-way Multidimensional Scaling analysis using Carroll and Chang (1970)'s INdividual Differences SCALing method as implemented using De Leeuw and Mair (2011)'s stress majorization algorithm.
indscal(..., type = "ordinal", itmax = 1000)
indscal(..., type = "ordinal", itmax = 1000)
... |
a sequence of datasets of class |
type |
is either "ratio" or "ordinal" |
itmax |
Maximum number of iterations |
an object of class INDSCAL
, i.e. a list containing
the following items:
delta
: Observed dissimilarities
obsdiss
: List of observed dissimilarities, normalized
confdiss
: List of configuration dissimilarities
conf
: List of matrices of final configurations
gspace
: Joint configurations aka group stimulus space
cweights
: Configuration weights
stress
: Stress-1 value
spp
: Stress per point
sps
: Stress per subject (matrix)
ndim
: Number of dimensions
model
: Type of smacof model
niter
: Number of iterations
nobj
: Number of objects
Jan de Leeuw and Patrick Mair
de Leeuw, J., & Mair, P. (2009). Multidimensional scaling using majorization: The R package smacof. Journal of Statistical Software, 31(3), 1-30, <https://www.jstatsoft.org/v31/i03/>
## Not run: attach(Namib) plot(indscal(DZ,HM,PT,Major,Trace)) ## End(Not run)
## Not run: attach(Namib) plot(indscal(DZ,HM,PT,Major,Trace)) ## End(Not run)
Turns a vector of numbers into an object of class KDE
using
a combination of the Botev (2010) bandwidth selector and the
Abramson (1982) adaptive kernel bandwidth modifier.
KDE(x, from = NA, to = NA, bw = NA, adaptive = TRUE, log = FALSE, n = 512, ...)
KDE(x, from = NA, to = NA, bw = NA, adaptive = TRUE, log = FALSE, n = 512, ...)
x |
a vector of numbers |
from |
minimum age of the time axis. If |
to |
maximum age of the time axis. If |
bw |
the bandwidth of the KDE. If NULL, |
adaptive |
boolean flag controlling if the adaptive KDE modifier of Abramson (1982) is used |
log |
transform the ages to a log scale if |
n |
horizontal resolution of the density estimate |
... |
optional arguments to be passed on to |
an object of class KDE
, i.e. a list containing the
following items:
x
: horizontal plot coordinates
y
: vertical plot coordinates
bw
: the base bandwidth of the density estimate
ages
: the data values from the input to the KDE
function
KDEs
data(Namib) samp <- Namib$DZ$x[['N1']] dens <- KDE(samp,0,3000,kernel="epanechnikov") plot(dens)
data(Namib) samp <- Namib$DZ$x[['N1']] dens <- KDE(samp,0,3000,kernel="epanechnikov") plot(dens)
KDEs
Convert a dataset of class distributional
into an object of
class KDEs
for further processing by the summaryplot
function.
KDEs( x, from = NA, to = NA, bw = NA, samebandwidth = TRUE, adaptive = TRUE, normalise = FALSE, log = FALSE, n = 512, ... )
KDEs( x, from = NA, to = NA, bw = NA, samebandwidth = TRUE, adaptive = TRUE, normalise = FALSE, log = FALSE, n = 512, ... )
x |
an object of class |
from |
minimum limit of the x-axis. |
to |
maximum limit of the x-axis. |
bw |
the bandwidth of the kernel density estimates. If
|
samebandwidth |
boolean flag indicating whether the same
bandwidth should be used for all samples. If
|
adaptive |
boolean flag switching on the adaptive bandwidth modifier of Abramson (1982) |
normalise |
boolean flag indicating whether or not the KDEs should all integrate to the same value. |
log |
boolean flag indicating whether the data should by plotted on a logarithmic scale. |
n |
horizontal resolution of the density estimates |
... |
optional parameters to be passed on to |
an object of class KDEs
, i.e. a list containing the
following items:
kdes
: a named list with objects of class KDE
from
: the beginning of the common time scale
to
: the end of the common time scale
themax
: the maximum probability density of all the KDEs
pch
: the plot symbol to be used by plot.KDEs
xlabel
: the x-axis label to be used by plot.KDEs
KDE
data(Namib) KDEs <- KDEs(Namib$DZ,0,3000,pch=NA) summaryplot(KDEs,ncol=3)
data(Namib) KDEs <- KDEs(Namib$DZ,0,3000,pch=NA) summaryplot(KDEs,ncol=3)
Returns the Kolmogorov-Smirnov dissimilarity between two samples
KS.diss(x, ...) ## Default S3 method: KS.diss(x, y, ...) ## S3 method for class 'distributional' KS.diss(x, ...)
KS.diss(x, ...) ## Default S3 method: KS.diss(x, y, ...) ## S3 method for class 'distributional' KS.diss(x, ...)
x |
the first sample as a vector |
... |
optional arguments |
y |
the second sample as a vector |
a scalar value representing the maximum vertical distance between the two cumulative distributions
data(Namib) print(KS.diss(Namib$DZ$x[['N1']],Namib$DZ$x[['T8']]))
data(Namib) print(KS.diss(Namib$DZ$x[['N1']],Namib$DZ$x[['T8']]))
Returns the Kuiper dissimilarity between two samples
Kuiper.diss(x, ...) ## Default S3 method: Kuiper.diss(x, y, ...) ## S3 method for class 'distributional' Kuiper.diss(x, ...)
Kuiper.diss(x, ...) ## Default S3 method: Kuiper.diss(x, y, ...) ## S3 method for class 'distributional' Kuiper.diss(x, ...)
x |
the first sample as a vector |
... |
optional arguments |
y |
the second sample as a vector |
a scalar value representing the sum of the maximum vertical distances above and below the cumulative distributions of x and y
data(Namib) print(Kuiper.diss(Namib$DZ$x[['N1']],Namib$DZ$x[['T8']]))
data(Namib) print(Kuiper.diss(Namib$DZ$x[['N1']],Namib$DZ$x[['T8']]))
Add lines to an existing ternary diagram
## S3 method for class 'ternary' lines(x, ...)
## S3 method for class 'ternary' lines(x, ...)
x |
an object of class |
... |
optional arguments to the generic |
tern <- ternary(Namib$PT,'Q',c('KF','P'),c('Lm','Lv','Ls')) plot(tern,pch=21,bg='red',labels=NULL) middle <- matrix(c(0.01,0.49,0.01,0.49,0.98,0.02),2,3) lines(ternary(middle))
tern <- ternary(Namib$PT,'Q',c('KF','P'),c('Lm','Lv','Ls')) plot(tern,pch=21,bg='red',labels=NULL) middle <- matrix(c(0.01,0.49,0.01,0.49,0.98,0.02),2,3) lines(ternary(middle))
Performs classical or nonmetric Multidimensional Scaling analysis of provenance data
MDS(x, ...) ## Default S3 method: MDS(x, classical = FALSE, k = 2, ...) ## S3 method for class 'compositional' MDS(x, classical = FALSE, k = 2, ...) ## S3 method for class 'counts' MDS(x, classical = FALSE, k = 2, ...) ## S3 method for class 'distributional' MDS(x, classical = FALSE, k = 2, nb = 0, ...) ## S3 method for class 'varietal' MDS(x, classical = FALSE, k = 2, nb = 0, ...)
MDS(x, ...) ## Default S3 method: MDS(x, classical = FALSE, k = 2, ...) ## S3 method for class 'compositional' MDS(x, classical = FALSE, k = 2, ...) ## S3 method for class 'counts' MDS(x, classical = FALSE, k = 2, ...) ## S3 method for class 'distributional' MDS(x, classical = FALSE, k = 2, nb = 0, ...) ## S3 method for class 'varietal' MDS(x, classical = FALSE, k = 2, nb = 0, ...)
x |
an object of class |
... |
optional arguments If If If If Otherwise, |
classical |
boolean flag indicating whether classical
( |
k |
the desired dimensionality of the solution |
nb |
number of bootstrap resamples. If |
an object of class MDS
, i.e. a list containing the
following items:
points
: a two column vector of the fitted configuration
classical
: a boolean flag indicating whether the MDS
configuration was obtained by classical (TRUE
) or nonmetric
(FALSE
) MDS.
diss
: the dissimilarity matrix used for the MDS analysis
stress
: (only if classical=TRUE
) the final stress
achieved (in percent)
Nordsvan, A.R., Kirscher, U., Kirkland, C.L., Barham, M. and Brennan, D.T., 2020. Resampling (detrital) zircon age distributions for accurate multidimensional scaling solutions. Earth-Science Reviews, p.103149.
Vermeesch, P., 2013, Multi-sample comparison of detrital age distributions. Chemical Geology v.341, 140-146, doi:10.1016/j.chemgeo.2013.01.010
data(Namib) plot(MDS(Namib$Major,classical=TRUE))
data(Namib) plot(MDS(Namib$Major,classical=TRUE))
Models grain size distribution of minerals and rock fragments of different densities
minsorting( X, dens, sname = NULL, phi = 2, sigmaphi = 1, medium = "freshwater", from = -2.25, to = 5.5, by = 0.25 )
minsorting( X, dens, sname = NULL, phi = 2, sigmaphi = 1, medium = "freshwater", from = -2.25, to = 5.5, by = 0.25 )
X |
an object of class |
dens |
a vector of mineral and rock densities |
sname |
sample name if unspecified, the first sample of the dataset will be used |
phi |
the mean grain size of the sample in Krumbein's phi units |
sigmaphi |
the standard deviation of the grain size distirbution, in phi units |
medium |
the transport medium, one of either "air", "freshwater" or "seawater" |
from |
the minimum grain size to be evaluated, in phi units |
to |
the maximum grain size to be evaluated, in phi units |
by |
the grain size interval of the output table, in phi units |
an object of class minsorting
, i.e. a list with two
tables:
mfract: the grain size distribution of each mineral (sum of the columns = 1)
mcomp: the composition of each grain size fraction (sum of the rows = 1)
Alberto Resentini and Pieter Vermeesch
Resentini, A, Malusa, M G and Garzanti, E. "MinSORTING: An Excel worksheet for modelling mineral grain-size distribution in sediments, with application to detrital geochronology and provenance studies." Computers & Geosciences 59 (2013): 90-97.
Garzanti, E, Ando, S and Vezzoli, G. "Settling equivalence of detrital minerals and grain-size dependence of sediment composition." Earth and Planetary Science Letters 273.1 (2008): 138-151.
restore
data(endmembers,densities) distribution <- minsorting(endmembers,densities,sname='ophiolite',phi=2, sigmaphi=1,medium="seawater",by=0.05) plot(distribution,cumulative=FALSE)
data(endmembers,densities) distribution <- minsorting(endmembers,densities,sname='ophiolite',phi=2, sigmaphi=1,medium="seawater",by=0.05) plot(distribution,cumulative=FALSE)
A large dataset of provenance data from Namibia comprised of 14 sand samples from the Namib Sand Sea and 2 samples from the Orange River.
Namib
is a list containing the following 6 items:
DZ
: a distributional
dataset containing the zircon
U-Pb ages for ca. 100 grains from each sample, as well as their
(1-sigma) analytical uncertainties.
PT
: a compositional
dataset with the bulk petrography
of the samples, i.e. the quartz (‘Q’), K-feldspar (‘KF’),
plagioclase (‘P’), and lithic fragments of metamorphic (‘Lm’),
volcanic (‘Lv’) and sedimentary (‘Ls’) origin.
HM
: a compositional
dataset containing the heavy
mineral composition of the samples, comprised of zircon (‘zr’),
tourmaline (‘tm’), rutile (‘rt’), Ti-oxides (‘TiOx’), titanite
(‘sph’), apatite (‘ap’), epidote (‘ep’), garnet (‘gt’), staurolite
(‘st’), andalusite (‘and’), kyanite (‘ky’), sillimanite (‘sil’),
amphibole (‘amp’), clinopyroxene (‘cpx’) and orthopyroxene (‘opx’).
PTHM
: a compositional
dataset combining the variables
contained in PT
and HM
plus ‘mica’, ‘opaques’,
'turbids' and 'other' transparent heavy minerals (‘LgM’),
normalised to 100.
Major
: a compositional
dataset listing the
concentrations (in wt
TiO2, P2O5 and MnO.
Trace
: a compositional
data listing the concentrations
(in ppm) of Rb, Sr, Ba, Sc, Y, La, Ce, Pr, Nd, Sm, Gd, Dy, Er, Yb, Th,
U, Zr, Hf, V, Nb, Cr, Co, Ni, Cu, Zn, Ga and Pb.
Pieter Vermeesch and Eduardo Garzanti
Vermeesch, P. and Garzanti, E., Making geological sense of 'Big Data' in sedimentary provenance analysis, Chemical Geology 409 (2015) 20-27
samp <- Namib$DZ$x[['N1']] dens <- KDE(samp,0,3000) plot(dens)
samp <- Namib$DZ$x[['N1']] dens <- KDE(samp,0,3000) plot(dens)
Performs PCA of compositional data using a centred logratio distance
PCA(x, ...)
PCA(x, ...)
x |
an object of class |
... |
optional arguments to R's |
an object of classes PCA
, which is synonymous to the
stats package's prcomp
class.
data(Namib) plot(MDS(Namib$Major,classical=TRUE)) dev.new() plot(PCA(Namib$Major),asp=1) print("This example demonstrates the equivalence of classical MDS and PCA")
data(Namib) plot(MDS(Namib$Major,classical=TRUE)) dev.new() plot(PCA(Namib$Major),asp=1) print("This example demonstrates the equivalence of classical MDS and PCA")
Plot the results of a correspondence analysis as a biplot
## S3 method for class 'CA' plot(x, labelcol = "black", vectorcol = "red", components = c(1, 2), ...)
## S3 method for class 'CA' plot(x, labelcol = "black", vectorcol = "red", components = c(1, 2), ...)
x |
an object of class |
labelcol |
colour of the sample labels (may be a vector). |
vectorcol |
colour of the vector loadings for the variables |
components |
two-element vector of components to be plotted |
... |
optional arguments of the generic |
CA
data(Namib) plot(CA(Namib$PT))
data(Namib) plot(CA(Namib$PT))
Plots an object of class compositional
as a pie chart
## S3 method for class 'compositional' plot(x, sname, annotate = TRUE, colmap = NULL, ...)
## S3 method for class 'compositional' plot(x, sname, annotate = TRUE, colmap = NULL, ...)
x |
an object of class |
sname |
the sample name |
annotate |
a boolean flag controlling if the pies of the pie-chart should be labelled |
colmap |
an optional string with the name of one of R's built-in colour palettes (e.g., heat.colors, terrain.colors, topo.colors, cm.colors), which are to be used for plotting the data. |
... |
optional parameters to be passed on to the graphics object |
data(Namib) plot(Namib$Major,'N1',colmap='heat.colors')
data(Namib) plot(Namib$Major,'N1',colmap='heat.colors')
Plot one or several samples from a distributional
dataset as
a histogram or Cumulative Age Distributions (CAD).
## S3 method for class 'distributional' plot( x, snames = NULL, annotate = TRUE, CAD = FALSE, pch = NA, verticals = TRUE, colmap = NULL, ... )
## S3 method for class 'distributional' plot( x, snames = NULL, annotate = TRUE, CAD = FALSE, pch = NA, verticals = TRUE, colmap = NULL, ... )
x |
an object of class |
snames |
a string or a vector of string with the names of the
samples that need plotting if |
annotate |
boolean flag indicating whether the x- and y-axis should be labelled |
CAD |
boolean flag indicating whether the data should be
plotted as a cumulative age distribution or a histogram. For
multi-sample plots, the function will override this value with
|
pch |
an optional symbol to mark the sample points along the CAD |
verticals |
boolean flag indicating if the horizontal lines of the CAD should be connected by vertical lines |
colmap |
an optional string with the name of one of R's built-in colour palettes (e.g., heat.colors, terrain.colors, topo.colors, cm.colors), which are to be used for plotting the data. |
... |
optional arguments to the generic |
data(Namib) plot(Namib$DZ,c('N1','N2'))
data(Namib) plot(Namib$DZ,c('N1','N2'))
Plots the group configuration of a Generalised Procrustes Analysis
## S3 method for class 'GPA' plot(x, pch = NA, pos = NULL, col = "black", bg = "white", cex = 1, ...)
## S3 method for class 'GPA' plot(x, pch = NA, pos = NULL, col = "black", bg = "white", cex = 1, ...)
x |
an object of class |
pch |
plot symbol |
pos |
position of the sample labels relative to the plot symbols if pch != NA |
col |
plot colour (may be a vector) |
bg |
background colour (may be a vector) |
cex |
relative size of plot symbols |
... |
optional arguments to the generic |
procrustes
data(Namib) GPA <- procrustes(Namib$DZ,Namib$HM) coast <- c('N1','N2','N3','N10','N11','N12','T8','T13') snames <- names(Namib$DZ) bgcol <- rep('yellow',length(snames)) bgcol[which(snames %in% coast)] <- 'red' plot(GPA,pch=21,bg=bgcol)
data(Namib) GPA <- procrustes(Namib$DZ,Namib$HM) coast <- c('N1','N2','N3','N10','N11','N12','T8','T13') snames <- names(Namib$DZ) bgcol <- rep('yellow',length(snames)) bgcol[which(snames %in% coast)] <- 'red' plot(GPA,pch=21,bg=bgcol)
Given an object of class INDSCAL
, generates two plots: the
group configuration and the subject weights. Together, these
describe a 3-way MDS model.
## S3 method for class 'INDSCAL' plot( x, asp = 1, pch = NA, pos = NULL, col = "black", bg = "white", cex = 1, xlab = "X", ylab = "Y", xaxt = "n", yaxt = "n", option = 2, ... )
## S3 method for class 'INDSCAL' plot( x, asp = 1, pch = NA, pos = NULL, col = "black", bg = "white", cex = 1, xlab = "X", ylab = "Y", xaxt = "n", yaxt = "n", option = 2, ... )
x |
an object of class |
asp |
the aspect ratio of the plot |
pch |
plot symbol (may be a vector) |
pos |
position of the sample labels relative to the plot symbols if pch != NA |
col |
plot colour (may be a vector) |
bg |
background colour (may be a vector) |
cex |
relative size of plot symbols |
xlab |
a string with the label of the x axis |
ylab |
a string with the label of the y axis |
xaxt |
if = 's', adds ticks to the x axis |
yaxt |
if = 's', adds ticks to the y axis |
option |
either:
|
... |
optional arguments to the generic plot function |
indscal
data(Namib) coast <- c('N1','N2','N3','N10','N11','N12','T8','T13') snames <- names(Namib$DZ) pch <- rep(21,length(snames)) pch[which(snames %in% coast)] <- 22 plot(indscal(Namib$DZ,Namib$HM),pch=pch)
data(Namib) coast <- c('N1','N2','N3','N10','N11','N12','T8','T13') snames <- names(Namib$DZ) pch <- rep(21,length(snames)) pch[which(snames %in% coast)] <- 22 plot(indscal(Namib$DZ,Namib$HM),pch=pch)
Plots an object of class KDE
## S3 method for class 'KDE' plot(x, pch = "|", xlab = "age [Ma]", ylab = "", ...)
## S3 method for class 'KDE' plot(x, pch = "|", xlab = "age [Ma]", ylab = "", ...)
x |
an object of class |
pch |
the symbol used to show the samples. May be a vector.
Set |
xlab |
the label of the x-axis |
ylab |
the label of the y-axis |
... |
optional parameters to be passed on to the graphics object |
KDE
data(Namib) samp <- Namib$DZ$x[['N1']] dens <- KDE(samp,from=0,to=3000) plot(dens)
data(Namib) samp <- Namib$DZ$x[['N1']] dens <- KDE(samp,from=0,to=3000) plot(dens)
Plots an object of class KDEs
## S3 method for class 'KDEs' plot(x, sname = NA, annotate = TRUE, pch = "|", ...)
## S3 method for class 'KDEs' plot(x, sname = NA, annotate = TRUE, pch = "|", ...)
x |
an object of class |
sname |
optional sample name. If |
annotate |
add a time axis? |
pch |
symbol to be used to mark the sample points along the
x-axis. Change to |
... |
optional parameters to be passed on to the
|
KDEs summaryplot
data(Namib) kdes <- KDEs(Namib$DZ) plot(kdes,ncol=2)
data(Namib) kdes <- KDEs(Namib$DZ) plot(kdes,ncol=2)
Plots the coordinates of a multidimensional scaling analysis as an X-Y scatter plot or ‘map’ and, if x$classical = FALSE, a Shepard plot.
## S3 method for class 'MDS' plot( x, nnlines = FALSE, pch = NA, pos = NULL, cex = 1, col = "black", bg = "white", oma = rep(1, 4), mar = rep(2, 4), mgp = c(2, 1, 0), xpd = NA, Shepard = 2, ... )
## S3 method for class 'MDS' plot( x, nnlines = FALSE, pch = NA, pos = NULL, cex = 1, col = "black", bg = "white", oma = rep(1, 4), mar = rep(2, 4), mgp = c(2, 1, 0), xpd = NA, Shepard = 2, ... )
x |
an object of class |
nnlines |
if TRUE, draws nearest neighbour lines |
pch |
plot character (see ?plot for details). May be a vector. |
pos |
position of the sample labels relative to the plot
symbols if |
cex |
relative size of plot symbols (see |
col |
plot colour (may be a vector) |
bg |
background colour (may be a vector) |
oma |
A vector of the form |
mar |
A numerical vector of the form |
mgp |
The margin line (in |
xpd |
A logical value or |
Shepard |
either:
|
... |
optional arguments to the generic |
MDS
data(Namib) mds <- MDS(Namib$DZ) coast <- c('N1','N2','N3','N10','N11','N12','T8','T13') snames <- names(Namib$DZ) bgcol <- rep('yellow',length(snames)) bgcol[which(snames %in% coast)] <- 'red' plot(mds,pch=21,bg=bgcol)
data(Namib) mds <- MDS(Namib$DZ) coast <- c('N1','N2','N3','N10','N11','N12','T8','T13') snames <- names(Namib$DZ) bgcol <- rep('yellow',length(snames)) bgcol[which(snames %in% coast)] <- 'red' plot(mds,pch=21,bg=bgcol)
Plot the grain size distributions of the different minerals under consideration
## S3 method for class 'minsorting' plot(x, cumulative = FALSE, components = NULL, ...)
## S3 method for class 'minsorting' plot(x, cumulative = FALSE, components = NULL, ...)
x |
an object of class |
cumulative |
boolean flag indicating whether the grain size distribution should be plotted as a density or cumulative probability curve. |
components |
string or list of strings with the names of a subcomposition that needs plotting |
... |
optional parameters to be passed on to graphics::matplot
(see |
minsorting
data(endmembers,densities) OPH <- subset(endmembers,select="ophiolite") distribution <- minsorting(OPH,densities,phi=2,sigmaphi=1, medium="air",by=0.05) plot(distribution,components=c('F','px','opaques'))
data(endmembers,densities) OPH <- subset(endmembers,select="ophiolite") distribution <- minsorting(OPH,densities,phi=2,sigmaphi=1, medium="air",by=0.05) plot(distribution,components=c('F','px','opaques'))
Plot the results of a principal components analysis as a biplot
## S3 method for class 'PCA' plot( x, labelcol = "black", vectorcol = "red", choices = 1L:2L, scale = 1, pc.biplot = FALSE, ... )
## S3 method for class 'PCA' plot( x, labelcol = "black", vectorcol = "red", choices = 1L:2L, scale = 1, pc.biplot = FALSE, ... )
x |
an object of class |
labelcol |
colour(s) of the sample labels (may be a vector). |
vectorcol |
colour of the vector loadings for the variables |
choices |
see the help pages of the generic |
scale |
see the help pages of the generic |
pc.biplot |
see the help pages of the generic |
... |
optional arguments of the generic |
PCA
data(Namib) plot(PCA(Namib$Major))
data(Namib) plot(PCA(Namib$Major))
Plots triplets of compositional data on a ternary diagram
## S3 method for class 'ternary' plot( x, type = "grid", pch = NA, pos = NULL, labels = names(x), showpath = FALSE, bg = NA, col = "cornflowerblue", ticks = seq(0, 1, 0.25), ticklength = 0.02, lty = 2, lwd = 1, ... )
## S3 method for class 'ternary' plot( x, type = "grid", pch = NA, pos = NULL, labels = names(x), showpath = FALSE, bg = NA, col = "cornflowerblue", ticks = seq(0, 1, 0.25), ticklength = 0.02, lty = 2, lwd = 1, ... )
x |
an object of class |
type |
adds annotations to the ternary diagram, one of either
|
pch |
plot character, see |
pos |
position of the sample labels relative to the plot symbols if pch != NA |
labels |
vector of strings to be added to the plot symbols |
showpath |
if |
bg |
background colour for the plot symbols (may be a vector) |
col |
colour to be used for the background lines (if applicable) |
ticks |
vector of tick values between 0 and 1 |
ticklength |
number between 0 and 1 to mark the length of the ticks |
lty |
line type for the annotations (see |
lwd |
line thickness for the annotations |
... |
optional arguments to the generic |
ternary
data(Namib) tern <- ternary(Namib$PT,'Q',c('KF','P'),c('Lm','Lv','Ls')) plot(tern,type='QFL.descriptive',pch=21,bg='red',labels=NULL)
data(Namib) tern <- ternary(Namib$PT,'Q',c('KF','P'),c('Lm','Lv','Ls')) plot(tern,type='QFL.descriptive',pch=21,bg='red',labels=NULL)
Add points to an existing ternary diagram
## S3 method for class 'ternary' points(x, ...)
## S3 method for class 'ternary' points(x, ...)
x |
an object of class |
... |
optional arguments to the generic |
tern <- ternary(Namib$PT,'Q',c('KF','P'),c('Lm','Lv','Ls')) plot(tern,pch=21,bg='red',labels=NULL) # add the geometric mean composition as a yellow square: gmean <- ternary(exp(colMeans(log(tern$x)))) points(gmean,pch=22,bg='yellow')
tern <- ternary(Namib$PT,'Q',c('KF','P'),c('Lm','Lv','Ls')) plot(tern,pch=21,bg='red',labels=NULL) # add the geometric mean composition as a yellow square: gmean <- ternary(exp(colMeans(log(tern$x)))) points(gmean,pch=22,bg='yellow')
Given a number of input datasets, this function performs an MDS
analysis on each of these and the feeds the resulting
configurations into the GPA()
function.
procrustes(...)
procrustes(...)
... |
a sequence of datasets of classes |
an object of class GPA
, i.e. a list containing the
following items:
points
: a two column vector with the coordinates of the
group configuration
labels
: a list with the sample names
Pieter Vermeesch
Gower, J.C. (1975). Generalized Procrustes analysis, Psychometrika, 40, 33-50.
GPA
data(Namib) gpa1 <- procrustes(Namib$DZ,Namib$HM) plot(gpa1) data(SNSM) gpa2 <- procrustes(SNSM$ap) plot(gpa2)
data(Namib) gpa1 <- procrustes(Namib$DZ,Namib$HM) plot(gpa1) data(SNSM) gpa2 <- procrustes(SNSM$ap) plot(gpa2)
provenance
For those less familiar with the syntax of the R
programming
language, the provenance()
function provides a user-friendly
way to access the most important functionality in the form of a
menu-based query interface. Further details and examples are
provided on https://www.ucl.ac.uk/~ucfbpve/provenance/
provenance
provides statistical tools to interpret large
amounts of distributional (single grain analyses) and compositional
(mineralogical and bulk chemical) data from the command line, or
using a menu-based user interface.
provenance()
provenance()
A list of documented functions may be viewed by typing
help(package='provenance')
. Detailed instructions are
provided at https://www.ucl.ac.uk/~ucfbpve/provenance/ and in
the Sedimentary Geology paper by Vermeesch, Resentini and Garzanti
(2016).
Pieter Vermeesch
Maintainer: Pieter Vermeesch [email protected]
Vermeesch, P., Resentini, A. and Garzanti, E., an R package for statistical provenance analysis, Sedimentary Geology, doi:10.1016/j.sedgeo.2016.01.009.
Vermeesch, P., Resentini, A. and Garzanti, E., 2016, An R package for statistical provenance analysis, Sedimentary Geology, 336, 14-25.
Useful links:
Implementation of a graphical device developed by Rex Galbraith to display several estimates of the same quantity that have different standard errors.
## S3 method for class 'counts' radialplot( x, num = 1, den = 2, from = NA, to = NA, t0 = NA, sigdig = 2, show.numbers = FALSE, pch = 21, levels = NA, clabel = "", bg = c("white", "red"), title = TRUE, ... )
## S3 method for class 'counts' radialplot( x, num = 1, den = 2, from = NA, to = NA, t0 = NA, sigdig = 2, show.numbers = FALSE, pch = 21, levels = NA, clabel = "", bg = c("white", "red"), title = TRUE, ... )
x |
an object of class |
num |
index or name of the numerator variable |
den |
index or name of the denominator variable |
from |
minimum limit of the radial scale |
to |
maximum limit of the radial scale |
t0 |
central value |
sigdig |
the number of significant digits of the numerical values reported in the title of the graphical output. |
show.numbers |
boolean flag ( |
pch |
plot character (default is a filled circle) |
levels |
a vector with additional values to be displayed as different background colours of the plot symbols. |
clabel |
label of the colour legend |
bg |
a vector of two background colours for the plot symbols.
If |
title |
add a title to the plot? |
... |
additional arguments to the generic |
The radial plot (Galbraith, 1988, 1990) is a graphical device that
was specifically designed to display heteroscedastic data, and is
constructed as follows. Consider a set of dates
and uncertainties
. Define
to be a transformation of
(e.g.,
),
and let
be its propagated analytical uncertainty
(i.e.,
in the case of a logarithmic
transformation). Create a scatterplot of
values,
where
and
,
where
is some reference value such as the mean. The
slope of a line connecting the origin of this scatterplot with any
of the
s is proportional to
and, hence,
the date
. These dates can be more easily visualised by
drawing a radial scale at some convenient distance from the origin
and annotating it with labelled ticks at the appropriate
angles. While the angular position of each data point represents
the date, its horizontal distance from the origin is proportional
to the precision. Imprecise measurements plot on the left hand side
of the radial plot, whereas precise age determinations are found
further towards the right. Thus, radial plots allow the observer to
assess both the magnitude and the precision of quantitative data in
one glance.
Galbraith, R.F., 1988. Graphical display of estimates having differing standard errors. Technometrics, 30(3), pp.271-281.
Galbraith, R.F., 1990. The radial plot: graphical assessment of spread in ages. International Journal of Radiation Applications and Instrumentation. Part D. Nuclear Tracks and Radiation Measurements, 17(3), pp.207-214.
Galbraith, R.F. and Laslett, G.M., 1993. Statistical models for mixed fission track ages. Nuclear Tracks and Radiation Measurements, 21(4), pp.459-470.
data(Namib) radialplot(Namib$PT,num='Q',den='P')
data(Namib) radialplot(Namib$PT,num='Q',den='P')
Reads a data table containing compositional data (e.g. chemical concentrations)
read.compositional( fname, method = NULL, colmap = "rainbow", sep = ",", dec = ".", row.names = 1, header = TRUE, check.names = FALSE, ... )
read.compositional( fname, method = NULL, colmap = "rainbow", sep = ",", dec = ".", row.names = 1, header = TRUE, check.names = FALSE, ... )
fname |
a string with the path to the .csv file |
method |
either "bray" (for the Bray-Curtis distance) or "aitchison" (for Aitchison's central logratio distance). If omitted, the function defaults to 'aitchison', unless there are zeros present in the data. |
colmap |
an optional string with the name of one of R's built-in colour palettes (e.g., heat.colors, terrain.colors, topo.colors, cm.colors), which are to be used for plotting the data. |
sep |
the field separator character. Values on each line of the file are separated by this character. |
dec |
the character used in the file for decimal points. |
row.names |
a vector of row names. This can be a vector giving the actual row names, or a single number giving the column of the which contains the row names, or character string the name of the table column containing the row names. |
header |
a logical value indicating whether the file contains the names of the variables as its first line. |
check.names |
logical. If |
... |
optional arguments to the built-in |
an object of class compositional
, i.e. a list with
the following items:
x
: a data frame with the samples as rows and the categories as columns
method
: either "aitchison" (for Aitchison's centred logratio
distance) or "bray" (for the Bray-Curtis distance)
colmap
: the colour map provided by the input argument
name
: the name of the data object, extracted from the file path
fname <- system.file("Namib/Major.csv",package="provenance") Major <- read.compositional(fname) plot(PCA(Major))
fname <- system.file("Namib/Major.csv",package="provenance") Major <- read.compositional(fname) plot(PCA(Major))
Reads a data table containing point-counting data (e.g. petrographic, heavy mineral, palaeontological or palynological data)
read.counts( fname, method = "chisq", colmap = "rainbow", sep = ",", dec = ".", row.names = 1, header = TRUE, check.names = FALSE, ... )
read.counts( fname, method = "chisq", colmap = "rainbow", sep = ",", dec = ".", row.names = 1, header = TRUE, check.names = FALSE, ... )
fname |
a string with the path to the .csv file |
method |
either "chisq" (for the chi-square distance) or "bray" (for the Bray-Curtis distance) |
colmap |
an optional string with the name of one of R's built-in colour palettes (e.g., heat.colors, terrain.colors, topo.colors, cm.colors), which are to be used for plotting the data. |
sep |
the field separator character. Values on each line of the file are separated by this character. |
dec |
the character used in the file for decimal points. |
row.names |
a vector of row names. This can be a vector giving the actual row names, or a single number giving the column of the which contains the row names, or character string the name of the table column containing the row names. |
header |
a logical value indicating whether the file contains the names of the variables as its first line. |
check.names |
logical. If |
... |
optional arguments to the built-in |
an object of class counts
, i.e. a list with the
following items:
x
: a data frame with the samples as rows and the categories
as columns
colmap
: the colour map provided by the input argument
name
: the name of the data object, extracted from the file path
fname <- system.file("Namib/HM.csv",package="provenance") Major <- read.counts(fname) #plot(PCA(HM))
fname <- system.file("Namib/HM.csv",package="provenance") Major <- read.counts(fname) #plot(PCA(HM))
Reads a data table containing densities to be used for hydraulic sorting corrections (minsorting and srd functions)
read.densities( fname, sep = ",", dec = ".", header = TRUE, check.names = FALSE, ... )
read.densities( fname, sep = ",", dec = ".", header = TRUE, check.names = FALSE, ... )
fname |
a string with the path to the .csv file |
sep |
the field separator character. Values on each line of the file are separated by this character. |
dec |
the character used in the file for decimal points. |
header |
a logical value indicating whether the file contains the names of the variables as its first line. |
check.names |
logical. If |
... |
optional arguments to the built-in |
a vector with mineral and rock densities
data(Namib,densities) N8 <- subset(Namib$HM,select="N8") distribution <- minsorting(N8,densities,phi=2,sigmaphi=1,medium="air",by=0.05) plot(distribution)
data(Namib,densities) N8 <- subset(Namib$HM,select="N8") distribution <- minsorting(N8,densities,phi=2,sigmaphi=1,medium="air",by=0.05) plot(distribution)
Reads a data table containing distributional data, i.e. lists of continuous data such as detrital zircon U-Pb ages.
read.distributional( fname, errorfile = NA, method = "KS", xlab = "age [Ma]", colmap = "rainbow", sep = ",", dec = ".", header = TRUE, check.names = FALSE, ... )
read.distributional( fname, errorfile = NA, method = "KS", xlab = "age [Ma]", colmap = "rainbow", sep = ",", dec = ".", header = TRUE, check.names = FALSE, ... )
fname |
the path of a .csv file with the input data, arranged in columns. |
errorfile |
the (optional) path of a .csv file with the
standard errors of the input data, arranged by column in the
same order as |
method |
an optional string specifying the dissimilarity
measure which should be used for comparing this with other
datasets. Should be one of either |
xlab |
an optional string specifying the nature and units of the data. This string is used to label kernel density estimates. |
colmap |
an optional string with the name of one of R's built-in colour palettes (e.g., heat.colors, terrain.colors, topo.colors, cm.colors), which are to be used for plotting the data. |
sep |
the field separator character. Values on each line of the file are separated by this character. |
dec |
the character used in the file for decimal points. |
header |
a logical value indicating whether the file contains the names of the variables as its first line. |
check.names |
logical. If |
... |
optional arguments to the built-in |
an object of class distributional
, i.e. a list with
the following items:
x
: a named list of vectors containing the numerical data for
each sample
err
: an (optional) named list of vectors containing the
standard errors of x
method
: either "KS" (for Kolmogorov-Smirnov), "Kuiper" (for
the Kuiper statistic) or "SH" (for Sircombe Hazelton)
breaks
: a vector with the locations of the histogram bin edges
xlab
: a string containing the label to be given to the
x-axis on all plots
colmap
: the colour map provided by the input argument
name
: the name of the data object, extracted from the file path
agefile <- system.file("Namib/DZ.csv",package="provenance") errfile <- system.file("Namib/DZerr.csv",package="provenance") DZ <- read.distributional(agefile,errfile) plot(KDE(DZ$x$N1))
agefile <- system.file("Namib/DZ.csv",package="provenance") errfile <- system.file("Namib/DZerr.csv",package="provenance") DZ <- read.distributional(agefile,errfile) plot(KDE(DZ$x$N1))
Reads a data table containing compositional data (e.g. chemical concentrations) for multiple grains and multiple samples
read.varietal( fname, snames = NULL, sep = ",", dec = ".", method = "KS", check.names = FALSE, row.names = 1, ... )
read.varietal( fname, snames = NULL, sep = ",", dec = ".", method = "KS", check.names = FALSE, row.names = 1, ... )
fname |
file name (character string) |
snames |
either a vector of sample names, an integer marking
the length of the sample name prefix, or
|
sep |
the field separator character. Values on each line of the file are separated by this character. |
dec |
the character used in the file for decimal points. |
method |
an optional string specifying the dissimilarity
measure which should be used for comparing this with other
datasets. Should be one of either |
check.names |
logical. If |
row.names |
logical. See the documentation for the
|
... |
optional arguments to the built-in |
an object of class varietal
, i.e. a list with the
following items:
x
: a compositional data table
snames
: a vector of strings corresponding to the sample names
name
: the name of the dataset, extracted from the file path
fn <- system.file("SNSM/Ttn_chem.csv",package="provenance") Ttn <- read.varietal(fname=fn,snames=3) plot(MDS(Ttn))
fn <- system.file("SNSM/Ttn_chem.csv",package="provenance") Ttn <- read.varietal(fname=fn,snames=3) plot(MDS(Ttn))
Restore the detrital composition back to a specified source rock density (SRD)
restore(X, dens, target = 2.71)
restore(X, dens, target = 2.71)
X |
an object of class |
dens |
a vector of rock and mineral densities |
target |
the target density (in g/cm3) |
an object of class SRDcorrected
, i.e. an object of class
compositional
which is a daughter of class compositional
containing the restored composition, plus one additional member called
restoration
, containing the intermediate steps of the SRD correction
algorithm.
Alberto Resentini and Pieter Vermeesch
Garzanti E, Ando, S and Vezzoli, G. "Settling equivalence of detrital minerals and grain-size dependence of sediment composition." Earth and Planetary Science Letters 273.1 (2008): 138-151.
minsorting
data(Namib,densities) rescomp <- restore(Namib$PTHM,densities,2.71) HMcomp <- c("zr","tm","rt","sph","ap","ep","gt", "st","amp","cpx","opx") amcomp <- amalgamate(rescomp,Plag="P",HM=HMcomp,Opq="opaques") plot(ternary(amcomp),showpath=TRUE)
data(Namib,densities) rescomp <- restore(Namib$PTHM,densities,2.71) HMcomp <- c("zr","tm","rt","sph","ap","ep","gt", "st","amp","cpx","opx") amcomp <- amalgamate(rescomp,Plag="P",HM=HMcomp,Opq="opaques") plot(ternary(amcomp),showpath=TRUE)
Calculates Sircombe and Hazelton's L2 distance between the Kernel Functional Estimates (KFEs, not to be confused with Kernel Density Estimates!) of two samples with specified analytical uncertainties
SH.diss(x, i, j, c.con = 0)
SH.diss(x, i, j, c.con = 0)
x |
an object of class |
i |
index of the first sample |
j |
index of the second sample |
c.con |
smoothing bandwidth of the kernel functional estimate |
a scalar value expressing the L2 distance between the KFEs of samples i and j
Keith Sircombe and Martin Hazelton
Sircombe, K. N., and M. L. Hazelton. "Comparison of detrital zircon age distributions by kernel functional estimation." Sedimentary Geology 171.1 (2004): 91-111.
KS.diss
datfile <- system.file("Namib/DZ.csv",package="provenance") errfile <- system.file("Namib/DZerr.csv",package="provenance") DZ <- read.distributional(datfile,errfile) d <- SH.diss(DZ,1,2) print(d)
datfile <- system.file("Namib/DZ.csv",package="provenance") errfile <- system.file("Namib/DZerr.csv",package="provenance") DZ <- read.distributional(datfile,errfile) d <- SH.diss(DZ,1,2) print(d)
A list of varietal datasets including detrital zircon
(zr
), apatite (ap
) and titanite
(tit
) compositions from the Sierra Nevada de Santa
Marta, provided by L. Caracciolo (FAU Erlangen).
Luca Caracciolo, Diana Hatzenbuehler and David Chew.
plot(MDS(SNSM$tit))
plot(MDS(SNSM$tit))
Return a subset of provenance data according to some specified indices
## S3 method for class 'distributional' subset(x, subset = NULL, select = NULL, ...) ## S3 method for class 'compositional' subset(x, subset = NULL, components = NULL, select = NULL, ...) ## S3 method for class 'counts' subset(x, subset = NULL, components = NULL, select = NULL, ...) ## S3 method for class 'varietal' subset(x, subset = NULL, components = NULL, select = NULL, ...)
## S3 method for class 'distributional' subset(x, subset = NULL, select = NULL, ...) ## S3 method for class 'compositional' subset(x, subset = NULL, components = NULL, select = NULL, ...) ## S3 method for class 'counts' subset(x, subset = NULL, components = NULL, select = NULL, ...) ## S3 method for class 'varietal' subset(x, subset = NULL, components = NULL, select = NULL, ...)
x |
an object of class |
subset |
logical expression indicating elements or rows to keep: missing values are taken as false. |
select |
a vector of sample names |
... |
optional arguments for the generic subset function |
components |
vector of categories (column names) to keep |
an object of the same class as x
data(Namib) coast <- c("N1","N2","T8","T13","N12","N13") ZTRcoast <- subset(Namib$HM,select=coast,components=c('gt','cpx','ep')) DZcoast <- subset(Namib$DZ,select=coast) summaryplot(ZTRcoast,KDEs(DZcoast),ncol=2)
data(Namib) coast <- c("N1","N2","T8","T13","N12","N13") ZTRcoast <- subset(Namib$HM,select=coast,components=c('gt','cpx','ep')) DZcoast <- subset(Namib$DZ,select=coast) summaryplot(ZTRcoast,KDEs(DZcoast),ncol=2)
Arranges kernel density estimates and pie charts in a grid format
summaryplot(..., ncol = 1, pch = NA)
summaryplot(..., ncol = 1, pch = NA)
... |
a sequence of datasets of class |
ncol |
the number of columns |
pch |
(optional) symbol to be used to mark the sample points along the x-axis of the KDEs (if appropriate). |
a summary plot of all the data comprised of KDEs for the
datasets of class KDEs
, pie charts for those of class
compositional
or counts
and histograms for those
of class distributional
.
data(Namib) KDEs <- KDEs(Namib$DZ,0,3000) summaryplot(KDEs,Namib$HM,Namib$PT,ncol=2)
data(Namib) KDEs <- KDEs(Namib$DZ,0,3000) summaryplot(KDEs,Namib$HM,Namib$PT,ncol=2)
Create an object of class ternary
ternary(X, x = 1, y = 2, z = 3)
ternary(X, x = 1, y = 2, z = 3)
X |
an object of class |
x |
string/number or a vector of strings/numbers indicating the variables/indices making up the first subcomposition of the ternary system. |
y |
second (set of) variables |
z |
third (set of) variables |
an object of class ternary
, i.e. a list containing:
x: a three column matrix (or vector) of ternary compositions.
and (if X is of class SRDcorrected
)
restoration: a list of intermediate ternary compositions inherited from the SRD correction
restore
data(Namib) tern <- ternary(Namib$PT,c('Q'),c('KF','P'),c('Lm','Lv','Ls')) plot(tern,type="QFL")
data(Namib) tern <- ternary(Namib$PT,c('Q'),c('KF','P'),c('Lm','Lv','Ls')) plot(tern,type="QFL")
plot a confidence region around the data or
around its mean.
ternary.ellipse(x, ...) ## Default S3 method: ternary.ellipse(x, alpha = 0.05, population = TRUE, ...) ## S3 method for class 'compositional' ternary.ellipse(x, alpha = 0.05, population = TRUE, ...) ## S3 method for class 'counts' ternary.ellipse(x, alpha = 0.05, population = TRUE, ...)
ternary.ellipse(x, ...) ## Default S3 method: ternary.ellipse(x, alpha = 0.05, population = TRUE, ...) ## S3 method for class 'compositional' ternary.ellipse(x, alpha = 0.05, population = TRUE, ...) ## S3 method for class 'counts' ternary.ellipse(x, alpha = 0.05, population = TRUE, ...)
x |
an object of class |
... |
optional formatting arguments |
alpha |
cutoff level for the confidence ellipse |
population |
show the standard deviation of the entire population or the standard error of the mean? |
data(Namib) tern <- ternary(Namib$Major,'CaO','Na2O','K2O') plot(tern) ternary.ellipse(tern)
data(Namib) tern <- ternary(Namib$Major,'CaO','Na2O','K2O') plot(tern) ternary.ellipse(tern)
Add text an existing ternary diagram
## S3 method for class 'ternary' text(x, labels = 1:nrow(x$x), ...)
## S3 method for class 'ternary' text(x, labels = 1:nrow(x$x), ...)
x |
an object of class |
labels |
a character vector or expression specifying the text to be written |
... |
optional arguments to the generic |
data(Namib) tern <- ternary(Namib$Major,'CaO','Na2O','K2O') plot(tern,pch=21,bg='red',labels=NULL) # add the geometric mean composition as a text label: gmean <- ternary(exp(colMeans(log(tern$x)))) text(gmean,labels='geometric mean')
data(Namib) tern <- ternary(Namib$Major,'CaO','Na2O','K2O') plot(tern,pch=21,bg='red',labels=NULL) # add the geometric mean composition as a text label: gmean <- ternary(exp(colMeans(log(tern$x)))) text(gmean,labels='geometric mean')
Convert an object of class varietal
either to a list of
distributional objects by breaking it up into separate elements, or
to a single distributional object corresponding to the first
principal component.
varietal2distributional(x, bycol = FALSE, plot = FALSE)
varietal2distributional(x, bycol = FALSE, plot = FALSE)
x |
an object of class |
bycol |
logical. If |
plot |
logical. If |
Ttn_file <- system.file("SNSM/Ttn_chem.csv",package="provenance") Ttn <- read.varietal(fn=Ttn_file,snames=3) varietal2distributional(Ttn,bycol=FALSE,plot=TRUE)
Ttn_file <- system.file("SNSM/Ttn_chem.csv",package="provenance") Ttn <- read.varietal(fn=Ttn_file,snames=3) varietal2distributional(Ttn,bycol=FALSE,plot=TRUE)
Returns the Wasserstein distance between two samples
Wasserstein.diss(x, ...) ## Default S3 method: Wasserstein.diss(x, y, ...) ## S3 method for class 'distributional' Wasserstein.diss(x, log = FALSE, ...) ## S3 method for class 'varietal' Wasserstein.diss(x, package = "transport", verbose = FALSE, ...)
Wasserstein.diss(x, ...) ## Default S3 method: Wasserstein.diss(x, y, ...) ## S3 method for class 'distributional' Wasserstein.diss(x, log = FALSE, ...) ## S3 method for class 'varietal' Wasserstein.diss(x, package = "transport", verbose = FALSE, ...)
x |
the first sample as a vector |
... |
optional arguments to the
|
y |
the second sample as a vector |
log |
logical. Take the lograthm of the data before calculating the distances? |
package |
the name of the package that provides the 2D
Wasserstein distance. Currently, this can be either
|
verbose |
logical. If |
a scalar value
The default S3 method was written by Pieter Vermeesch,
using modified code from Dominic Schuhmacher's transport
package (transport1d
function), as implemented in
IsoplotR
.
data(Namib) print(Wasserstein.diss(Namib$DZ$x[['N1']],Namib$DZ$x[['T8']]))
data(Namib) print(Wasserstein.diss(Namib$DZ$x[['N1']],Namib$DZ$x[['T8']]))