Title: | Multivariate Functional Principal Component Analysis for Data Observed on Different Dimensional Domains |
---|---|
Description: | Calculate a multivariate functional principal component analysis for data observed on different dimensional domains. The estimation algorithm relies on univariate basis expansions for each element of the multivariate functional data (Happ & Greven, 2018) <doi:10.1080/01621459.2016.1273115>. Multivariate and univariate functional data objects are represented by S4 classes for this type of data implemented in the package 'funData'. For more details on the general concepts of both packages and a case study, see Happ-Kurz (2020) <doi:10.18637/jss.v093.i05>. |
Authors: | Clara Happ-Kurz [aut, cre] |
Maintainer: | Clara Happ-Kurz <[email protected]> |
License: | GPL-2 |
Version: | 1.3-9 |
Built: | 2024-10-26 04:41:34 UTC |
Source: | https://github.com/clarahapp/mfpca |
This function implements the functional CP-TPA (FCP-TPA) algorithm, that
calculates a smooth PCA for 3D tensor data (i.e. N
observations of 2D
images with dimension S1 x S2
). The results are given in a
CANDECOMP/PARAFRAC (CP) model format
where
stands for the outer product,
is a scalar and
are eigenvectors for each direction of the tensor. In
this representation, the outer product
can
be regarded as the
-th eigenimage, while
represents the vector of individual scores for this eigenimage and each
observation.
FCP_TPA( X, K, penMat, alphaRange, verbose = FALSE, tol = 1e-04, maxIter = 15, adaptTol = TRUE )
FCP_TPA( X, K, penMat, alphaRange, verbose = FALSE, tol = 1e-04, maxIter = 15, adaptTol = TRUE )
X |
The data tensor of dimensions |
K |
The number of eigentensors to be calculated. |
penMat |
A list with entries |
alphaRange |
A list of length 2 with entries |
verbose |
Logical. If |
tol |
A numeric value, giving the tolerance for relative error values in
the algorithm. Defaults to |
maxIter |
A numeric value, the maximal iteration steps. Can be doubled,
if |
adaptTol |
Logical. If |
The smoothness of the eigenvectors is induced by penalty
matrices for both image directions, that are weighted by smoothing parameters
. The eigenvectors
are not smoothed,
hence the algorithm does not induce smoothness along observations.
Optimal smoothing parameters are found via a nested generalized cross
validation. In each iteration of the TPA (tensor power algorithm), the GCV
criterion is optimized via optimize
on the interval
specified via alphaRange$v
(or alphaRange$w
, respectively).
The FCP_TPA algorithm is an iterative algorithm. Convergence is assumed if
the relative difference between the actual and the previous values are all
below the tolerance level tol
. The tolerance level is increased
automatically, if the algorithm has not converged after maxIter
steps
and if adaptTol = TRUE
. If the algorithm did not converge after
maxIter
steps (or 2 * maxIter
) steps, the function throws a
warning.
d |
A vector of length |
U |
A matrix of dimensions |
V |
A
matrix of dimensions |
W |
A matrix of dimensions |
G. I. Allen, "Multi-way Functional Principal Components Analysis", IEEE International Workshop on Computational Advances in Multi-Sensor Adaptive Processing, 2013.
# set.seed(1234) N <- 100 S1 <- 75 S2 <- 75 # define "true" components v <- sin(seq(-pi, pi, length.out = S1)) w <- exp(seq(-0.5, 1, length.out = S2)) # simulate tensor data with dimensions N x S1 x S2 X <- rnorm(N, sd = 0.5) %o% v %o% w # create penalty matrices (penalize first differences for each dimension) Pv <- crossprod(diff(diag(S1))) Pw <- crossprod(diff(diag(S2))) # estimate one eigentensor res <- FCP_TPA(X, K = 1, penMat = list(v = Pv, w = Pw), alphaRange = list(v = c(1e-4, 1e4), w = c(1e-4, 1e4)), verbose = TRUE) # plot the results and compare to true values plot(res$V) points(v/sqrt(sum(v^2)), pch = 20) legend("topleft", legend = c("True", "Estimated"), pch = c(20, 1)) plot(res$W) points(w/sqrt(sum(w^2)), pch = 20) legend("topleft", legend = c("True", "Estimated"), pch = c(20, 1))
# set.seed(1234) N <- 100 S1 <- 75 S2 <- 75 # define "true" components v <- sin(seq(-pi, pi, length.out = S1)) w <- exp(seq(-0.5, 1, length.out = S2)) # simulate tensor data with dimensions N x S1 x S2 X <- rnorm(N, sd = 0.5) %o% v %o% w # create penalty matrices (penalize first differences for each dimension) Pv <- crossprod(diff(diag(S1))) Pw <- crossprod(diff(diag(S2))) # estimate one eigentensor res <- FCP_TPA(X, K = 1, penMat = list(v = Pv, w = Pw), alphaRange = list(v = c(1e-4, 1e4), w = c(1e-4, 1e4)), verbose = TRUE) # plot the results and compare to true values plot(res$V) points(v/sqrt(sum(v^2)), pch = 20) legend("topleft", legend = c("True", "Estimated"), pch = c(20, 1)) plot(res$W) points(w/sqrt(sum(w^2)), pch = 20) legend("topleft", legend = c("True", "Estimated"), pch = c(20, 1))
This function calculates a multivariate functional principal component
analysis (MFPCA) based on i.i.d. observations of a
multivariate functional data-generating process
with elements
defined on a domain
. In particular, the
elements can be defined on different (dimensional) domains. The results
contain the mean function, the estimated multivariate functional principal
components
(having the same structure
as
), the associated eigenvalues
and the individual scores
. Moreover,
estimated trajectories for each observation based on the truncated
Karhunen-Loeve representation
are given
if desired (fit = TRUE
). The implementation of the observations
, the mean function and
multivariate functional principal components
uses the
multiFunData
class, which is defined
in the package funData.
MFPCA( mFData, M, uniExpansions, weights = rep(1, length(mFData)), fit = FALSE, approx.eigen = FALSE, bootstrap = FALSE, nBootstrap = NULL, bootstrapAlpha = 0.05, bootstrapStrat = NULL, verbose = options()$verbose )
MFPCA( mFData, M, uniExpansions, weights = rep(1, length(mFData)), fit = FALSE, approx.eigen = FALSE, bootstrap = FALSE, nBootstrap = NULL, bootstrapAlpha = 0.05, bootstrapStrat = NULL, verbose = options()$verbose )
mFData |
A |
M |
The number of multivariate functional principal components to calculate. |
uniExpansions |
A list characterizing the (univariate) expansion that is calculated for each element. See Details. |
weights |
An optional vector of weights, defaults to |
fit |
Logical. If |
approx.eigen |
Logical. If |
bootstrap |
Logical. If |
nBootstrap |
The number of bootstrap iterations to use. Defaults to
|
bootstrapAlpha |
A vector of numerics (or a single number) giving the
significance level for bootstrap intervals. Defaults to |
bootstrapStrat |
A stratification variable for bootstrap. Must be a
factor of length |
verbose |
Logical. If |
If the elements vary considerably in domain,
range or variation, a weight vector can be supplied
and the MFPCA is based on the weighted scalar product
and the
corresponding weighted covariance operator .
If bootstrap = TRUE
, pointwise bootstrap
confidence bands are generated for the multivariate eigenvalues as well as for multivariate functional principal
components
. The parameter
nBootstrap
gives the number of bootstrap
iterations. In each iteration, the observations are resampled on the level of
(multivariate) functions and the whole MFPCA is recalculated. In particular,
if the univariate basis depends on the data (FPCA approaches), basis
functions and scores are both re-estimated. If the basis functions are fixed
(e.g. splines), the scores from the original estimate are used to speed up
the calculations. The confidence bands for the eigenfunctions are calculated
separately for each element as pointwise percentile bootstrap confidence
intervals. Analogously, the confidence bands for the eigenvalues are also
percentile bootstrap confidence bands. The significance level(s) can be
defined by the bootstrapAlpha
parameter, which defaults to 5%. As a
result, the MFPCA
function returns a list CI
of the same length
as bootstrapAlpha
, containing the lower and upper bounds of the
confidence bands for the principal components as multiFunData
objects
of the same structure as mFData
. The confidence bands for the
eigenvalues are returned in a list CIvalues
, containing the upper and
lower bounds for each significance level.
The multivariate functional principal
component analysis relies on a univariate basis expansion for each element
. The univariate basis representation is calculated using
the
univDecomp
function, that passes the univariate functional
observations and optional parameters to the specific function. The univariate
decompositions are specified via the uniExpansions
argument in the
MFPCA
function. It is a list of the same length as the mFData
object, i.e. having one entry for each element of the multivariate functional
data. For each element, uniExpansion
must specify at least the type of
basis functions to use. Additionally, one may add further parameters. The
following basis representations are supported:
Given basis
functions. Then uniExpansions[[j]] = list(type = "given", functions,
scores, ortho)
, where functions
is a funData
object on the
same domain as mFData
, containing the given basis functions. The
parameters scores
and ortho
are optional. scores
is an
N x K
matrix containing the scores (or coefficients) of the observed
functions for the given basis functions, where N
is the number of
observed functions and K
is the number of basis functions. Note that
the scores need to be demeaned to give meaningful results. If scores are not
supplied, they are calculated using the given basis functions. The parameter
ortho
specifies whether the given basis functions are orthonormal
orhto = TRUE
or not ortho = FALSE
. If ortho
is not
supplied, the functions are treated as non-orthogonal. scores
and
ortho
are not checked for plausibility, use them at your own risk!
Univariate functional principal component analysis. Then
uniExpansions[[j]] = list(type = "uFPCA", nbasis, pve, npc, makePD)
,
where nbasis,pve,npc,makePD
are parameters passed to the
PACE
function for calculating the univariate functional
principal component analysis.
Basis functions expansions from the
package fda. Then uniExpansions[[j]] = list(type = "fda", ...)
,
where ...
are passed to funData2fd
, which
heavily builds on eval.fd
. If fda is not available,
a warning is thrown.
Spline basis functions (not penalized). Then
uniExpansions[[j]] = list(type = "splines1D", bs, m, k)
, where
bs,m,k
are passed to the functions univDecomp
and
univExpansion
. For two-dimensional tensor product splines, use
type = "splines2D"
.
Spline basis functions (with smoothness
penalty). Then uniExpansions[[j]] = list(type = "splines1Dpen", bs, m,
k)
, where bs,m,k
are passed to the functions univDecomp
and univExpansion
. Analogously to the unpenalized case, use
type = "splines2Dpen"
for 2D penalized tensor product splines.
Cosine basis functions. Use uniExpansions[[j]] = list(type = "DCT2D",
qThresh, parallel)
for functions one two-dimensional domains (images) and
type = "DCT3D"
for 3D images. The calculation is based on the discrete
cosine transform (DCT) implemented in the C-library fftw3
. If this
library is not available, the function will throw a warning. qThresh
gives the quantile for hard thresholding the basis coefficients based on
their absolute value. If parallel = TRUE
, the coefficients for
different images are calculated in parallel.
See univDecomp
and univExpansion
for details.
An object of class MFPCAfit
containing the following
components:
values |
A vector of estimated eigenvalues |
functions |
A
|
scores |
A matrix of dimension |
vectors |
A matrix representing the eigenvectors associated with the combined univariate score vectors. This might be helpful for calculating predictions. |
normFactors |
The normalizing factors used for calculating the multivariate eigenfunctions and scores. This might be helpful when calculation predictions. |
meanFunction |
A multivariate functional
data object, corresponding to the mean function. The MFPCA is applied to
the de-meaned functions in |
fit |
A
|
CI |
A
list of the same length as |
CIvalues |
A list of the same
length as |
C. Happ, S. Greven (2018): Multivariate Functional Principal Component Analysis for Data Observed on Different (Dimensional) Domains. Journal of the American Statistical Association, 113(522): 649-659. DOI: doi:10.1080/01621459.2016.1273115
C. Happ-Kurz (2020): Object-Oriented Software for Functional Data. Journal of Statistical Software, 93(5): 1-38. DOI: doi:10.18637/jss.v093.i05
See Happ-Kurz (2020. doi:10.18637/jss.v093.i05) for a general
introduction to the funData package and it's interplay with
MFPCA. This file also includes a case study on how to use
MFPCA
. Useful functions: multiFunData
,
PACE
, univDecomp
, univExpansion
,
summary
,
plot
,
scoreplot
oldPar <- par(no.readonly = TRUE) set.seed(1) ### simulate data (one-dimensional domains) sim <- simMultiFunData(type = "split", argvals = list(seq(0,1,0.01), seq(-0.5,0.5,0.02)), M = 5, eFunType = "Poly", eValType = "linear", N = 100) # MFPCA based on univariate FPCA uFPCA <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "uFPCA"), list(type = "uFPCA"))) summary(uFPCA) plot(uFPCA) # plot the eigenfunctions as perturbations of the mean scoreplot(uFPCA) # plot the scores # MFPCA based on univariate spline expansions splines <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "splines1D", k = 10), list(type = "splines1D", k = 10)), fit = TRUE) # calculate reconstruction, too summary(splines) plot(splines) # plot the eigenfunctions as perturbations of the mean scoreplot(splines) # plot the scores ### Compare estimates to true eigenfunctions # flip to make results more clear uFPCA$functions <- flipFuns(sim$trueFuns, uFPCA$functions) splines$functions <- flipFuns(sim$trueFuns, splines$functions) par(mfrow = c(1,2)) plot(sim$trueFuns[[1]], main = "Eigenfunctions\n1st Element", lwd = 2) plot(uFPCA$functions[[1]], lty = 2, add = TRUE) plot(splines$functions[[1]], lty = 3, add = TRUE) plot(sim$trueFuns[[2]], main = "Eigenfunctions\n2nd Element", lwd = 2) plot(uFPCA$functions[[2]], lty = 2, add = TRUE) plot(splines$functions[[2]], lty = 3, add = TRUE) legend("bottomleft", c("True", "uFPCA", "splines"), lty = 1:3, lwd = c(2,1,1)) # Test reconstruction for the first 10 observations plot(sim$simData[[1]], obs = 1:10, main = "Reconstruction\n1st Element", lwd = 2) plot(splines$fit[[1]], obs = 1:10, lty = 2, col = 1, add = TRUE) plot(sim$simData[[2]], obs = 1:10, main = "Reconstruction\n2nd Element", lwd = 2) plot(splines$fit[[2]], obs = 1:10, lty = 2, col = 1, add = TRUE) legend("bottomleft", c("True", "Reconstruction"), lty = c(1,2), lwd = c(2,1)) # MFPCA with Bootstrap-CI for the first 2 eigenfunctions ### ATTENTION: Takes long splinesBoot <- MFPCA(sim$simData, M = 2, uniExpansions = list(list(type = "splines1D", k = 10), list(type = "splines1D", k = 10)), bootstrap = TRUE, nBootstrap = 100, bootstrapAlpha = c(0.05, 0.1), verbose = TRUE) summary(splinesBoot) plot(splinesBoot$functions[[1]], ylim = c(-2,1.5)) plot(splinesBoot$CI$alpha_0.05$lower[[1]], lty = 2, add = TRUE) plot(splinesBoot$CI$alpha_0.05$upper[[1]], lty = 2, add = TRUE) plot(splinesBoot$CI$alpha_0.1$lower[[1]], lty = 3, add = TRUE) plot(splinesBoot$CI$alpha_0.1$upper[[1]], lty = 3, add = TRUE) abline(h = 0, col = "gray") plot(splinesBoot$functions[[2]], ylim = c(-1,2.5)) plot(splinesBoot$CI$alpha_0.05$lower[[2]], lty = 2, add = TRUE) plot(splinesBoot$CI$alpha_0.05$upper[[2]], lty = 2, add = TRUE) plot(splinesBoot$CI$alpha_0.1$lower[[2]], lty = 3, add = TRUE) plot(splinesBoot$CI$alpha_0.1$upper[[2]], lty = 3, add = TRUE) abline(h = 0, col = "gray") legend("topleft", c("Estimate", "95% CI", "90% CI"), lty = 1:3, lwd = c(2,1,1)) # Plot 95% confidence bands for eigenvalues plot(1:2, splinesBoot$values, pch = 20, ylim = c(0, 1.5), main = "Estimated eigenvalues with 95% CI", xlab = "Eigenvalue no.", ylab = "") arrows(1:2, splinesBoot$CIvalues$alpha_0.05$lower, 1:2, splinesBoot$CIvalues$alpha_0.05$upper, length = 0.05, angle = 90, code = 3) points(1:2, sim$trueVals[1:2], pch = 20, col = 4) legend("topright", c("Estimate", "True value"), pch = 20, col = c(1,4)) ### simulate data (two- and one-dimensional domains) ### ATTENTION: Takes long set.seed(2) sim <- simMultiFunData(type = "weighted", argvals = list(list(seq(0,1,0.01), seq(-1,1,0.02)), list(seq(-0.5,0.5,0.01))), M = list(c(4,5), 20), eFunType = list(c("Fourier", "Fourier"), "Poly"), eValType = "exponential", N = 150) # MFPCA based on univariate spline expansions (for images) and univariate FPCA (for functions) pca <- MFPCA(sim$simData, M = 10, uniExpansions = list(list(type = "splines2D", k = c(10,12)), list(type = "uFPCA"))) summary(pca) plot(pca) # plot the eigenfunctions as perturbations of the mean scoreplot(pca) # plot the scores ### Compare to true eigenfunctions # flip to make results more clear pca$functions <- flipFuns(sim$trueFuns[1:10], pca$functions) par(mfrow = c(5,2), mar = rep(2,4)) for(m in 2:6) # for m = 1, image.plot (used in plot(funData)) produces an error... { plot(sim$trueFuns[[1]], main = paste("True, m = ", m), obs = m) plot(pca$functions[[1]], main = paste("Estimate, m = ", m), obs = m) } par(mfrow = c(1,1)) plot(sim$trueFuns[[2]], main = "Eigenfunctions (2nd element)", lwd = 2, obs= 1:5) plot(pca$functions[[2]], lty = 2, add = TRUE, obs= 1:5) legend("bottomleft", c("True", "MFPCA"), lty = 1:2, lwd = c(2,1)) par(oldPar)
oldPar <- par(no.readonly = TRUE) set.seed(1) ### simulate data (one-dimensional domains) sim <- simMultiFunData(type = "split", argvals = list(seq(0,1,0.01), seq(-0.5,0.5,0.02)), M = 5, eFunType = "Poly", eValType = "linear", N = 100) # MFPCA based on univariate FPCA uFPCA <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "uFPCA"), list(type = "uFPCA"))) summary(uFPCA) plot(uFPCA) # plot the eigenfunctions as perturbations of the mean scoreplot(uFPCA) # plot the scores # MFPCA based on univariate spline expansions splines <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "splines1D", k = 10), list(type = "splines1D", k = 10)), fit = TRUE) # calculate reconstruction, too summary(splines) plot(splines) # plot the eigenfunctions as perturbations of the mean scoreplot(splines) # plot the scores ### Compare estimates to true eigenfunctions # flip to make results more clear uFPCA$functions <- flipFuns(sim$trueFuns, uFPCA$functions) splines$functions <- flipFuns(sim$trueFuns, splines$functions) par(mfrow = c(1,2)) plot(sim$trueFuns[[1]], main = "Eigenfunctions\n1st Element", lwd = 2) plot(uFPCA$functions[[1]], lty = 2, add = TRUE) plot(splines$functions[[1]], lty = 3, add = TRUE) plot(sim$trueFuns[[2]], main = "Eigenfunctions\n2nd Element", lwd = 2) plot(uFPCA$functions[[2]], lty = 2, add = TRUE) plot(splines$functions[[2]], lty = 3, add = TRUE) legend("bottomleft", c("True", "uFPCA", "splines"), lty = 1:3, lwd = c(2,1,1)) # Test reconstruction for the first 10 observations plot(sim$simData[[1]], obs = 1:10, main = "Reconstruction\n1st Element", lwd = 2) plot(splines$fit[[1]], obs = 1:10, lty = 2, col = 1, add = TRUE) plot(sim$simData[[2]], obs = 1:10, main = "Reconstruction\n2nd Element", lwd = 2) plot(splines$fit[[2]], obs = 1:10, lty = 2, col = 1, add = TRUE) legend("bottomleft", c("True", "Reconstruction"), lty = c(1,2), lwd = c(2,1)) # MFPCA with Bootstrap-CI for the first 2 eigenfunctions ### ATTENTION: Takes long splinesBoot <- MFPCA(sim$simData, M = 2, uniExpansions = list(list(type = "splines1D", k = 10), list(type = "splines1D", k = 10)), bootstrap = TRUE, nBootstrap = 100, bootstrapAlpha = c(0.05, 0.1), verbose = TRUE) summary(splinesBoot) plot(splinesBoot$functions[[1]], ylim = c(-2,1.5)) plot(splinesBoot$CI$alpha_0.05$lower[[1]], lty = 2, add = TRUE) plot(splinesBoot$CI$alpha_0.05$upper[[1]], lty = 2, add = TRUE) plot(splinesBoot$CI$alpha_0.1$lower[[1]], lty = 3, add = TRUE) plot(splinesBoot$CI$alpha_0.1$upper[[1]], lty = 3, add = TRUE) abline(h = 0, col = "gray") plot(splinesBoot$functions[[2]], ylim = c(-1,2.5)) plot(splinesBoot$CI$alpha_0.05$lower[[2]], lty = 2, add = TRUE) plot(splinesBoot$CI$alpha_0.05$upper[[2]], lty = 2, add = TRUE) plot(splinesBoot$CI$alpha_0.1$lower[[2]], lty = 3, add = TRUE) plot(splinesBoot$CI$alpha_0.1$upper[[2]], lty = 3, add = TRUE) abline(h = 0, col = "gray") legend("topleft", c("Estimate", "95% CI", "90% CI"), lty = 1:3, lwd = c(2,1,1)) # Plot 95% confidence bands for eigenvalues plot(1:2, splinesBoot$values, pch = 20, ylim = c(0, 1.5), main = "Estimated eigenvalues with 95% CI", xlab = "Eigenvalue no.", ylab = "") arrows(1:2, splinesBoot$CIvalues$alpha_0.05$lower, 1:2, splinesBoot$CIvalues$alpha_0.05$upper, length = 0.05, angle = 90, code = 3) points(1:2, sim$trueVals[1:2], pch = 20, col = 4) legend("topright", c("Estimate", "True value"), pch = 20, col = c(1,4)) ### simulate data (two- and one-dimensional domains) ### ATTENTION: Takes long set.seed(2) sim <- simMultiFunData(type = "weighted", argvals = list(list(seq(0,1,0.01), seq(-1,1,0.02)), list(seq(-0.5,0.5,0.01))), M = list(c(4,5), 20), eFunType = list(c("Fourier", "Fourier"), "Poly"), eValType = "exponential", N = 150) # MFPCA based on univariate spline expansions (for images) and univariate FPCA (for functions) pca <- MFPCA(sim$simData, M = 10, uniExpansions = list(list(type = "splines2D", k = c(10,12)), list(type = "uFPCA"))) summary(pca) plot(pca) # plot the eigenfunctions as perturbations of the mean scoreplot(pca) # plot the scores ### Compare to true eigenfunctions # flip to make results more clear pca$functions <- flipFuns(sim$trueFuns[1:10], pca$functions) par(mfrow = c(5,2), mar = rep(2,4)) for(m in 2:6) # for m = 1, image.plot (used in plot(funData)) produces an error... { plot(sim$trueFuns[[1]], main = paste("True, m = ", m), obs = m) plot(pca$functions[[1]], main = paste("Estimate, m = ", m), obs = m) } par(mfrow = c(1,1)) plot(sim$trueFuns[[2]], main = "Eigenfunctions (2nd element)", lwd = 2, obs= 1:5) plot(pca$functions[[2]], lty = 2, add = TRUE, obs= 1:5) legend("bottomleft", c("True", "MFPCA"), lty = 1:2, lwd = c(2,1)) par(oldPar)
Calculate multivariate basis expansion
multivExpansion(multiFuns, scores)
multivExpansion(multiFuns, scores)
multiFuns |
A multivariate functional data object, containing the multivariate basis functions |
scores |
A matrix containing the scores for each observation in each row. The number of columns must match the number of basis functions. |
A multiFunData
object containing the expanded functions for
each observation.
This function calculates a univariate functional principal components
analysis by smoothed covariance based on code from
fpca.sc
in package refund.
PACE( funDataObject, predData = NULL, nbasis = 10, pve = 0.99, npc = NULL, makePD = FALSE, cov.weight.type = "none" )
PACE( funDataObject, predData = NULL, nbasis = 10, pve = 0.99, npc = NULL, makePD = FALSE, cov.weight.type = "none" )
funDataObject |
An object of class |
predData |
An object of class |
nbasis |
An integer, representing the number of B-spline basis
functions used for estimation of the mean function and bivariate smoothing
of the covariance surface. Defaults to |
pve |
A numeric value between 0 and 1, the proportion of variance
explained: used to choose the number of principal components. Defaults to
|
npc |
An integer, giving a prespecified value for the number of
principal components. Defaults to |
makePD |
Logical: should positive definiteness be enforced for the
covariance surface estimate? Defaults to |
cov.weight.type |
The type of weighting used for the smooth covariance
estimate. Defaults to |
mu |
A |
values |
A vector containing the estimated eigenvalues. |
functions |
A
|
scores |
An matrix of estimated scores for the
observations in |
fit |
A |
npc |
The number of functional
principal components: either the supplied |
sigma2 |
The estimated measurement error variance (cf.
|
estVar |
The estimated smooth variance function of the data. |
This function works only for univariate functional data observed on one-dimensional domains.
funData
,
fpcaBasis
, univDecomp
oldPar <- par(no.readonly = TRUE) # simulate data sim <- simFunData(argvals = seq(-1,1,0.01), M = 5, eFunType = "Poly", eValType = "exponential", N = 100) # calculate univariate FPCA pca <- PACE(sim$simData, npc = 5) # Plot the results par(mfrow = c(1,2)) plot(sim$trueFuns, lwd = 2, main = "Eigenfunctions") # flip estimated functions for correct signs plot(flipFuns(sim$trueFuns,pca$functions), lty = 2, add = TRUE) legend("bottomright", c("True", "Estimate"), lwd = c(2,1), lty = c(1,2)) plot(sim$simData, lwd = 2, main = "Some Observations", obs = 1:7) plot(pca$fit, lty = 2, obs = 1:7, add = TRUE) # estimates are almost equal to true values legend("bottomright", c("True", "Estimate"), lwd = c(2,1), lty = c(1,2)) par(oldPar)
oldPar <- par(no.readonly = TRUE) # simulate data sim <- simFunData(argvals = seq(-1,1,0.01), M = 5, eFunType = "Poly", eValType = "exponential", N = 100) # calculate univariate FPCA pca <- PACE(sim$simData, npc = 5) # Plot the results par(mfrow = c(1,2)) plot(sim$trueFuns, lwd = 2, main = "Eigenfunctions") # flip estimated functions for correct signs plot(flipFuns(sim$trueFuns,pca$functions), lty = 2, add = TRUE) legend("bottomright", c("True", "Estimate"), lwd = c(2,1), lty = c(1,2)) plot(sim$simData, lwd = 2, main = "Some Observations", obs = 1:7) plot(pca$fit, lty = 2, obs = 1:7, add = TRUE) # estimates are almost equal to true values legend("bottomright", c("True", "Estimate"), lwd = c(2,1), lty = c(1,2)) par(oldPar)
Plots the eigenfunctions as perturbations of the mean (i.e. the mean function plus/minus a constant factor times each eigenfunction separately). If all elements have a one-dimensional domain, the plots can be combined, otherwise the effects of adding and subtracting are shown in two separate rows for each eigenfunction.
## S3 method for class 'MFPCAfit' plot( x, plotPCs = seq_len(nObs(x$functions)), stretchFactor = NULL, combined = FALSE, ... )
## S3 method for class 'MFPCAfit' plot( x, plotPCs = seq_len(nObs(x$functions)), stretchFactor = NULL, combined = FALSE, ... )
x |
An object of class |
plotPCs |
The principal components to be plotted. Defaults to all
components in the |
stretchFactor |
The factor by which the principal components are
multiplied before adding / subtracting them from the mean function. If
|
combined |
Logical: Should the plots be combined? (Works only if all
dimensions are one-dimensional). Defaults to |
... |
Further graphical parameters passed to the plot.funData functions for functional data. |
A plot of the principal components as perturbations of the mean.
# Simulate multivariate functional data on one-dimensonal domains # and calculate MFPCA (cf. MFPCA help) set.seed(1) # simulate data (one-dimensional domains) sim <- simMultiFunData(type = "split", argvals = list(seq(0,1,0.01), seq(-0.5,0.5,0.02)), M = 5, eFunType = "Poly", eValType = "linear", N = 100) # MFPCA based on univariate FPCA PCA <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "uFPCA"), list(type = "uFPCA"))) # Plot the results plot(PCA, combined = TRUE) # combine addition and subtraction in one plot
# Simulate multivariate functional data on one-dimensonal domains # and calculate MFPCA (cf. MFPCA help) set.seed(1) # simulate data (one-dimensional domains) sim <- simMultiFunData(type = "split", argvals = list(seq(0,1,0.01), seq(-0.5,0.5,0.02)), M = 5, eFunType = "Poly", eValType = "linear", N = 100) # MFPCA based on univariate FPCA PCA <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "uFPCA"), list(type = "uFPCA"))) # Plot the results plot(PCA, combined = TRUE) # combine addition and subtraction in one plot
Predict functions based on a truncated multivariate Karhunen-Loeve representation:
with estimated mean function and principal components
. The scores
can be either estimated (reconstruction
of observed functions) or user-defined (construction of new functions).
## S3 method for class 'MFPCAfit' predict(object, scores = object$scores, ...)
## S3 method for class 'MFPCAfit' predict(object, scores = object$scores, ...)
object |
An object of class |
scores |
A matrix containing the score values. The number of columns in
|
... |
Arguments passed to or from other methods. |
A multiFunData
object containing the predicted functions.
#' # Simulate multivariate functional data on one-dimensonal domains # and calculate MFPCA (cf. MFPCA help) set.seed(1) # simulate data (one-dimensional domains) sim <- simMultiFunData(type = "split", argvals = list(seq(0,1,0.01), seq(-0.5,0.5,0.02)), M = 5, eFunType = "Poly", eValType = "linear", N = 100) # MFPCA based on univariate FPCA PCA <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "uFPCA"), list(type = "uFPCA"))) # Reconstruct the original data pred <- predict(PCA) # default reconstructs data used for the MFPCA fit # plot the results: 1st element plot(sim$simData[[1]]) # original data plot(pred[[1]], add = TRUE, lty = 2) # reconstruction # plot the results: 2nd element plot(sim$simData[[2]]) # original data plot(pred[[2]], add = TRUE, lty = 2) # reconstruction
#' # Simulate multivariate functional data on one-dimensonal domains # and calculate MFPCA (cf. MFPCA help) set.seed(1) # simulate data (one-dimensional domains) sim <- simMultiFunData(type = "split", argvals = list(seq(0,1,0.01), seq(-0.5,0.5,0.02)), M = 5, eFunType = "Poly", eValType = "linear", N = 100) # MFPCA based on univariate FPCA PCA <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "uFPCA"), list(type = "uFPCA"))) # Reconstruct the original data pred <- predict(PCA) # default reconstructs data used for the MFPCA fit # plot the results: 1st element plot(sim$simData[[1]]) # original data plot(pred[[1]], add = TRUE, lty = 2) # reconstruction # plot the results: 2nd element plot(sim$simData[[2]]) # original data plot(pred[[2]], add = TRUE, lty = 2) # reconstruction
A print
function for class MFPCAfit
.
## S3 method for class 'MFPCAfit' print(x, ...)
## S3 method for class 'MFPCAfit' print(x, ...)
x |
An object of class |
... |
Arguments passed to or from other methods. |
No return value, called for side effects
A print
method for class MFPCAfit.summary
## S3 method for class 'summary.MFPCAfit' print(x, ...)
## S3 method for class 'summary.MFPCAfit' print(x, ...)
x |
An object of class |
... |
Arguments passed to or from other methods. |
No return value, called for side effects
Redirects to plot.default
scoreplot(PCAobject, ...)
scoreplot(PCAobject, ...)
PCAobject |
A principal component object |
... |
Arguments passed from or to other methods |
A bivariate plot of scores.
This function plots two scores of a multivariate functional principal component analysis for each observation.
## S3 method for class 'MFPCAfit' scoreplot(PCAobject, choices = 1:2, scale = FALSE, ...)
## S3 method for class 'MFPCAfit' scoreplot(PCAobject, choices = 1:2, scale = FALSE, ...)
PCAobject |
An object of class |
choices |
The indices of the scores that should by displayed. Defaults
to |
scale |
Logical. Should the scores be scaled by the estimated
eigenvalues to emphasize the proportions of total variance explained by the
components. Defaults to |
... |
Further parameters passed to the
|
A bivariate plot of scores.
# and calculate MFPCA (cf. MFPCA help) set.seed(1) # simulate data (one-dimensional domains) sim <- simMultiFunData(type = "split", argvals = list(seq(0,1,0.01), seq(-0.5,0.5,0.02)), M = 5, eFunType = "Poly", eValType = "linear", N = 100) # MFPCA based on univariate FPCA PCA <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "uFPCA"), list(type = "uFPCA"))) # Plot the first two scores scoreplot(PCA) # no scaling (default) scoreplot(PCA, scale = TRUE) # scale the scores by the first two eigenvalues
# and calculate MFPCA (cf. MFPCA help) set.seed(1) # simulate data (one-dimensional domains) sim <- simMultiFunData(type = "split", argvals = list(seq(0,1,0.01), seq(-0.5,0.5,0.02)), M = 5, eFunType = "Poly", eValType = "linear", N = 100) # MFPCA based on univariate FPCA PCA <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "uFPCA"), list(type = "uFPCA"))) # Plot the first two scores scoreplot(PCA) # no scaling (default) scoreplot(PCA, scale = TRUE) # scale the scores by the first two eigenvalues
This function plots the proportion of variance explained by the leading eigenvalues in an MFPCA against the number of the principal component.
## S3 method for class 'MFPCAfit' screeplot( x, npcs = min(10, length(x$values)), type = "lines", ylim = NULL, main = deparse(substitute(x)), ... )
## S3 method for class 'MFPCAfit' screeplot( x, npcs = min(10, length(x$values)), type = "lines", ylim = NULL, main = deparse(substitute(x)), ... )
x |
An object of class MFPCAfit, typically returned by a call to
|
npcs |
The number of eigenvalued to be plotted. Defaults to all eigenvalues if their number is less or equal to 10, otherwise show only the leading first 10 eigenvalues. |
type |
The type of screeplot to be plotted. Can be either
|
ylim |
The limits for the y axis. Can be passed either as a vector
of length 2 or as |
main |
The title of the plot. Defaults to the variable name of
|
... |
Other graphic parameters passed to
|
A screeplot, showing the decrease of the principal component score.
# Simulate multivariate functional data on one-dimensonal domains # and calculate MFPCA (cf. MFPCA help) set.seed(1) # simulate data (one-dimensional domains) sim <- simMultiFunData(type = "split", argvals = list(seq(0,1,0.01), seq(-0.5,0.5,0.02)), M = 5, eFunType = "Poly", eValType = "linear", N = 100) # MFPCA based on univariate FPCA PCA <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "uFPCA"), list(type = "uFPCA"))) # screeplot screeplot(PCA) # default options screeplot(PCA, npcs = 3, type = "barplot", main= "Screeplot")
# Simulate multivariate functional data on one-dimensonal domains # and calculate MFPCA (cf. MFPCA help) set.seed(1) # simulate data (one-dimensional domains) sim <- simMultiFunData(type = "split", argvals = list(seq(0,1,0.01), seq(-0.5,0.5,0.02)), M = 5, eFunType = "Poly", eValType = "linear", N = 100) # MFPCA based on univariate FPCA PCA <- MFPCA(sim$simData, M = 5, uniExpansions = list(list(type = "uFPCA"), list(type = "uFPCA"))) # screeplot screeplot(PCA) # default options screeplot(PCA, npcs = 3, type = "barplot", main= "Screeplot")
A summary
method for class MFPCAfit
## S3 method for class 'MFPCAfit' summary(object, ...)
## S3 method for class 'MFPCAfit' summary(object, ...)
object |
An object of class |
... |
Arguments passed to or from other methods. |
An object of class summary.MFPCAfit
Functionality adapted from the MATLAB tensor toolbox (https://www.tensortoolbox.org/).
ttv(A, v, dim)
ttv(A, v, dim)
A |
An array. |
v |
A list of the same length as |
dim |
A vector specifying the dimensions for the multiplication. |
Let A
be a tensor with dimensions and let
v
be a vector of length
. Then the tensor-vector-product along the
-th dimension is
defined as
It can hence be seen as a generalization of the matrix-vector product.
The tensor-vector-product along several dimensions between a tensor A
and multiple vectors v_1,...,v_k
() is defined as a
series of consecutive tensor-vector-product along the different dimensions.
For consistency, the multiplications are calculated from the dimension of the
highest order to the lowest.
An array, the result of the multiplication.
B. W. Bader and T. G. Kolda. Algorithm 862: MATLAB tensor classes for fast algorithm prototyping, ACM Transactions on Mathematical Software 32(4):635-653, December 2006.
# create a three-mode tensor a1 <- seq(0,1, length.out = 10) a2 <- seq(-1,1, length.out = 20) a3 <- seq(-pi, pi, length.out = 15) A <-a1 %o% a2 %o% a3 dim(A) # multiply along different dimensions dim(ttv(A = A, v = list(rnorm(10)), dim = 1)) dim(ttv(A = A, v = list(rnorm(20)), dim = 2)) dim(ttv(A = A, v = list(rnorm(15)), dim = 3)) # multiply along more than one dimension length(ttv(A = A, v = list(rnorm(10), rnorm(15)), dim = c(1,3)))
# create a three-mode tensor a1 <- seq(0,1, length.out = 10) a2 <- seq(-1,1, length.out = 20) a3 <- seq(-pi, pi, length.out = 15) A <-a1 %o% a2 %o% a3 dim(A) # multiply along different dimensions dim(ttv(A = A, v = list(rnorm(10)), dim = 1)) dim(ttv(A = A, v = list(rnorm(20)), dim = 2)) dim(ttv(A = A, v = list(rnorm(15)), dim = 3)) # multiply along more than one dimension length(ttv(A = A, v = list(rnorm(10), rnorm(15)), dim = c(1,3)))
This function implements the uncorrelated multilinear principal component analysis for tensors of dimension 2, 3 or 4. The code is basically the same as in the MATLAB toolbox UMPCA by Haiping Lu (Link: https://www.mathworks.com/matlabcentral/fileexchange/35432-uncorrelated-multilinear-principal-component-analysis-umpca, see also references).
UMPCA(TX, numP)
UMPCA(TX, numP)
TX |
The input training data in tensorial representation, the last mode
is the sample mode. For |
numP |
The dimension of the projected vector, denoted as |
Us |
The multilinear projection, consisting of |
TXmean |
The mean
of the input training samples |
odrIdx |
The ordering index of projected features in decreasing variance. |
As this algorithm aims more at uncorrelated features than at an optimal reconstruction of the data, hence it might give poor results when used for the univariate decomposition of images in MFPCA.
Haiping Lu, K.N. Plataniotis, and A.N. Venetsanopoulos, "Uncorrelated Multilinear Principal Component Analysis for Unsupervised Multilinear Subspace Learning", IEEE Transactions on Neural Networks, Vol. 20, No. 11, Page: 1820-1836, Nov. 2009.
set.seed(12345) # define "true" components a <- sin(seq(-pi, pi, length.out = 100)) b <- exp(seq(-0.5, 1, length.out = 150)) # simulate tensor data X <- a %o% b %o% rnorm(80, sd = 0.5) # estimate one component UMPCAres <- UMPCA(X, numP = 1) # plot the results and compare to true values plot(UMPCAres$Us[[1]][,1]) points(a/sqrt(sum(a^2)), pch = 20) # eigenvectors are defined only up to a sign change! legend("topright", legend = c("True", "Estimated"), pch = c(20, 1)) plot(UMPCAres$Us[[2]][,1]) points(b/sqrt(sum(b^2)), pch = 20) legend("topleft", legend = c("True", "Estimated"), pch = c(20, 1))
set.seed(12345) # define "true" components a <- sin(seq(-pi, pi, length.out = 100)) b <- exp(seq(-0.5, 1, length.out = 150)) # simulate tensor data X <- a %o% b %o% rnorm(80, sd = 0.5) # estimate one component UMPCAres <- UMPCA(X, numP = 1) # plot the results and compare to true values plot(UMPCAres$Us[[1]][,1]) points(a/sqrt(sum(a^2)), pch = 20) # eigenvectors are defined only up to a sign change! legend("topright", legend = c("True", "Estimated"), pch = c(20, 1)) plot(UMPCAres$Us[[2]][,1]) points(b/sqrt(sum(b^2)), pch = 20) legend("topleft", legend = c("True", "Estimated"), pch = c(20, 1))
This function calculates a univariate basis decomposition for a (univariate) functional data object.
univDecomp(type, funDataObject, ...)
univDecomp(type, funDataObject, ...)
type |
A character string, specifying the basis for which the decomposition is to be calculated. |
funDataObject |
A |
... |
Further parameters, passed to the function for the particular basis to use. |
Functional data can often be approximated by a linear
combination of basis functions
The basis functions may be
prespecified (such as spline basis functions or Fourier bases) or can
be estimated from the data (e.g. by functional principal component
analysis) and are the same for all observations . The coefficients (or scores)
reflect the
weight of each basis function
for the observed function
and can be used to characterize the individual
observations.
scores |
A matrix of scores (coefficients) for each observation based on the prespecified basis functions. |
B |
A
matrix containing the scalar products of the basis functions. Can be
|
ortho |
Logical. If |
functions |
A functional data object, representing
the basis functions. Can be |
The options type = "DCT2D"
and type =
"DCT3D"
have not been tested with ATLAS/MKL/OpenBLAS.
MFPCA
, univExpansion
,
fpcaBasis
, splineBasis1D
,
splineBasis1Dpen
, splineBasis2D
,
splineBasis2Dpen
, umpcaBasis
,
fcptpaBasis
, fdaBasis
,
dctBasis2D
, dctBasis3D
# generate some data dat <- simFunData(argvals = seq(0,1,0.01), M = 5, eFunType = "Poly", eValType = "linear", N = 100)$simData # decompose the data in univariate functional principal components... decFPCA <- univDecomp(type = "uFPCA", funDataObject = dat, npc = 5) str(decFPCA) # or in splines (penalized) decSplines <- univDecomp(type = "splines1Dpen", funDataObject = dat) # use mgcv's default params str(decSplines)
# generate some data dat <- simFunData(argvals = seq(0,1,0.01), M = 5, eFunType = "Poly", eValType = "linear", N = 100)$simData # decompose the data in univariate functional principal components... decFPCA <- univDecomp(type = "uFPCA", funDataObject = dat, npc = 5) str(decFPCA) # or in splines (penalized) decSplines <- univDecomp(type = "splines1Dpen", funDataObject = dat) # use mgcv's default params str(decSplines)
This function calculates a univariate basis expansion based on given scores (coefficients) and basis functions.
univExpansion( type, scores, argvals = ifelse(!is.null(functions), functions@argvals, NULL), functions, params = NULL )
univExpansion( type, scores, argvals = ifelse(!is.null(functions), functions@argvals, NULL), functions, params = NULL )
type |
A character string, specifying the basis for which the decomposition is to be calculated. |
scores |
A matrix of scores (coefficients) for each observation based on the given basis functions. |
argvals |
A list, representing the domain of the basis functions.
If |
functions |
A functional data object, representing the basis
functions. Can be |
params |
A list containing the parameters for the particular basis to use. |
This function calculates functional data
that is represented as a linear combination of basis functions
The basis functions may be prespecified (such as spline
basis functions or Fourier bases) or can be estimated from observed
data (e.g. by functional principal component analysis). If type =
"default"
(i.e. a linear combination of arbitrary basis functions is
to be calculated), both scores and basis functions must be supplied.
An object of class funData
with N
observations on
argvals
, corresponding to the linear combination of the basis
functions.
The options type = "spline2Dpen"
, type =
"DCT2D"
and type = "DCT3D"
have not been tested with
ATLAS/MKL/OpenBLAS.
MFPCA
, splineFunction1D
,
splineFunction2D
, splineFunction2Dpen
,
dctFunction2D
, dctFunction3D
,
expandBasisFunction
oldPar <- par(no.readonly = TRUE) par(mfrow = c(1,1)) set.seed(1234) ### Spline basis ### # simulate coefficients (scores) for N = 10 observations and K = 8 basis functions N <- 10 K <- 8 scores <- t(replicate(n = N, rnorm(K, sd = (K:1)/K))) dim(scores) # expand spline basis on [0,1] funs <- univExpansion(type = "splines1D", scores = scores, argvals = list(seq(0,1,0.01)), functions = NULL, # spline functions are known, need not be given params = list(bs = "ps", m = 2, k = K)) # params for mgcv plot(funs, main = "Spline reconstruction") ### PCA basis ### # simulate coefficients (scores) for N = 10 observations and K = 8 basis functions N <- 10 K <- 8 scores <- t(replicate(n = N, rnorm(K, sd = (K:1)/K))) dim(scores) # Fourier basis functions as eigenfunctions eFuns <- eFun(argvals = seq(0,1,0.01), M = K, type = "Fourier") # expand eigenfunction basis funs <- univExpansion(type = "uFPCA", scores = scores, argvals = NULL, # use argvals of eFuns (default) functions = eFuns) plot(funs, main = "PCA reconstruction") par(oldPar)
oldPar <- par(no.readonly = TRUE) par(mfrow = c(1,1)) set.seed(1234) ### Spline basis ### # simulate coefficients (scores) for N = 10 observations and K = 8 basis functions N <- 10 K <- 8 scores <- t(replicate(n = N, rnorm(K, sd = (K:1)/K))) dim(scores) # expand spline basis on [0,1] funs <- univExpansion(type = "splines1D", scores = scores, argvals = list(seq(0,1,0.01)), functions = NULL, # spline functions are known, need not be given params = list(bs = "ps", m = 2, k = K)) # params for mgcv plot(funs, main = "Spline reconstruction") ### PCA basis ### # simulate coefficients (scores) for N = 10 observations and K = 8 basis functions N <- 10 K <- 8 scores <- t(replicate(n = N, rnorm(K, sd = (K:1)/K))) dim(scores) # Fourier basis functions as eigenfunctions eFuns <- eFun(argvals = seq(0,1,0.01), M = K, type = "Fourier") # expand eigenfunction basis funs <- univExpansion(type = "uFPCA", scores = scores, argvals = NULL, # use argvals of eFuns (default) functions = eFuns) plot(funs, main = "PCA reconstruction") par(oldPar)