| albatross-package {albatross} | R Documentation |
PARAFAC Analysis of Fluorescence Excitation-Emission Matrices
Description
Day after day, day after day,
We stuck, nor breath nor motion;
As idle as a painted ship
Upon a painted ocean.
Water, water, every where,
And all the boards did shrink;
Water, water, every where,
Nor any drop to drink.
– Samuel Taylor Coleridge, The Rime of the Ancient Mariner
Perform parallel factor analysis (PARAFAC: Hitchcock, 1927) <doi:10.1002/sapm192761164> on fluorescence excitation-emission matrices: handle scattering signal and inner filter effect, scale the dataset, fit the model; perform split-half validation or jack-knifing. Modified approaches such as Whittaker interpolation, randomised split-half, and fluorescence and scattering model estimation are also available. The package has a low dependency footprint and has been tested on a wide range of R versions.
Details
In order to work with your data, create feem and/or
feemcube objects from files or matrix or array objects.
Use feemlist to import files in bulk. If your files
aren't in one of the formats supported by feem but you
can read their contents by other means, you can supply an importer
function to feemlist; it should take a file name and
return the corresponding feem object.
Operations that can be performed on the objects include plotting
(plot.feem), calculation of fluorescence indices
(feemindex), inner-filter effect correction
(feemife), handling of scattering signal
(feemscatter), changing the wavelength grid of the data
by means of interpolation (feemgrid), and scaling
(feemscale). Scaling may be automatically undone after
performing the PARAFAC decomposition so that the resulting scores
would correspond to the data as it was before the scaling.
All processing functions can take individual feem
objects, lists of them, or feemcube objects and return
values of the appropriate kind. For example, feemscatter
always returns an object of the same class but with the scattering
signal handled, while feemindex returns named numeric
vectors for individual feems but
data.frames for collections of them. There's a
slight memory benefit to using lists of feem objects,
but the difference shouldn't be noticeable, so there's nothing to
worry about if you started with a feemcube.
In order to compute PARAFAC, you need to convert your data into a
feemcube. Whether you perform jack-knifing, split-half
analysis, or PARAFAC itself, a copy of the data cube is kept together
with the results and can be extracted back using the
feemcube function. The resulting objects support a
plot method (described in the same help page) and can give you
the data as a few-column data.frame using the
coef method.
Once the analysis is finished, the PARAFAC model can be exported for
the OpenFluor database (write.openfluor) or stored as an
R object using standard R tools (save or
saveRDS).
Index of help topics:
[.feem Extract or replace parts of FEEM objects
[.feemcube Extract or replace parts of FEEM cubes
absindex Functions of absorbance data
albatross-package PARAFAC Analysis of Fluorescence
Excitation-Emission Matrices
as.data.frame.feem Transform a FEEM object into a data.frame
feem Create a fluorescence excitation-emission
matrix object
feemcorcondia Core consistency diagnostic for PARAFAC models
feemcube Data cubes of fluorescence excitation-emission
matrices
feemflame Fluorescence and scAttering Model Estimation
feemgrid Interpolate FEEMs on a given wavelength grid
feemife Absorbance-based inner filter effect correction
feemindex Fluorescence indices and peak values
feemjackknife Jack-knife outlier detection in PARAFAC models
feemlist Create lists of FEEM objects
feemparafac Compute PARAFAC on a FEEM cube object and
access the results
feems Synthetic fluorescence excitation-emission
matrices and absorbance spectra
feemscale Rescale FEEM spectra to a given norm and
remember the scale factor
feemscatter Handle scattering signal in FEEMs
feemsplithalf Split-half analysis of PARAFAC models
marine.colours Perceptually uniform palettes
plot.feem Plot a FEEM object
write.openfluor Export a PARAFAC model for the OpenFluor
database
Author(s)
Ivan Krylov [aut, cre], Timur Labutin [ths], Anastasia Drozdova [rev]
References
Murphy KR, Stedmon CA, Graeber D, Bro R (2013). “Fluorescence spectroscopy and multi-way techniques. PARAFAC.” Analytical Methods, 5, 6557-6566. doi:10.1039/c3ay41160e.
Pucher M, Wünsch U, Weigelhofer G, Murphy K, Hein T, Graeber D (2019). “staRdom: Versatile Software for Analyzing Spectroscopic Data of Dissolved Organic Matter in R.” Water, 11(11), 2366. doi:10.3390/w11112366.
Cleese J, Jones T (1970). “Albatross: Flavours of different sea birds.” Journal of Flying Circus, 1.13, 7:05-7:45.
Krylov I, Drozdova A, Labutin T (2020). “Albatross R package to study PARAFAC components of DOM fluorescence from mixing zones of arctic shelf seas.” Chemometrics and Intelligent Laboratory Systems, 207(104176). doi:10.1016/j.chemolab.2020.104176.
See Also
feem, feemlist, feemindex,
feemife, feemscatter,
feemgrid, feemcube,
feemscale, feemsplithalf,
feemparafac, feemjackknife,
feemflame, absindex.
Examples
data(feems)
dataset <- feemcube(feems, FALSE)
dataset <- feemscatter(dataset, c(24, 15, 10), 'whittaker')
dataset <- feemife(dataset, absorp)
plot(dataset <- feemscale(dataset, na.rm = TRUE))
# takes a long time
(sh <- feemsplithalf(
feemscatter(cube, c(24, 15)),
# in real life, set a stricter stopping criterion, 1e-6..1e-8
nfac = 2:4, splits = 4, ctol = 1e-4
))
plot(sh)
pf <- feemparafac(cube, nfac = 3, ctol = 1e-4)
plot(pf)