Reference API
Documentation for Multispati.jl
's public interface.
Index
Multispati.AbstractMultispati
Multispati.MULTISPATI
Multispati.SpatialPCA
Base.size
LinearAlgebra.eigvals
LinearAlgebra.eigvecs
Multispati.moransIbounds
Multispati.varianceMoransIdecomposition
MultivariateStats.projection
MultivariateStats.reconstruct
MultivariateStats.reconstruct
Statistics.mean
StatsAPI.fit
StatsAPI.fit
StatsAPI.predict
StatsAPI.predict
API
Multispati
Multispati.AbstractMultispati
— TypeAbstractMultispati
Multispati.MULTISPATI
— TypeMULTISPATI <: AbstractMultispati
StatsAPI.fit
— Methodfit(MULTISPATI, X, W, Q=I, D=I / size(X, 2); ...)
Perform MULTISPATI over the data given a matrix X
. Each column of X
is an observation. W
is a connectivity matrix where $w_{ij}$ is the connection from j -> i. Q
is a symmetric matrix of size n
(or LinearAlgebra.UniformScaling(@ref)) and D
a symmetric matrix of size d
(or LinearAlgebra.UniformScaling(@ref))
Keyword arguments
maxoutdim
: The output dimension, i.e. dimension of the transformed space (min(d, nc-1))solver
: The choice of solver::eig
: usesLinearAlgebra.eigen
(default):eigs
: usesArpack.eigs
(always used for sparse data)
tol
: Convergence tolerance foreigs
solver (default0.0
)maxiter
: Maximum number of iterations foreigs
solver (default300
)
References
StatsAPI.predict
— Methodpredict(M::MULTISPATI, x::AbstractVecOrMat{<:Real})
Transform the observations x
with the model M
.
Here, x
can be either a vector of length d
or a matrix where each column is an observation.
MultivariateStats.reconstruct
— Methodreconstruct(M::MULTISPATI, y::AbstractVecOrMat{<:Real})
Approximately reconstruct the observations y
to the original space using the model M
.
Here, y
can be either a vector of length p
or a matrix where each column gives the components for an observation.
Multispati.moransIbounds
— FunctionmoransIbounds(M::AbstractMultispati; sparse_approx::Bool=true)
Return the bounds and expected value for Moran's I given the model M
in the order $I_{min}$, $I_{max}$, $I_0$
Base.size
— Methodsize(M)
Returns a tuple with the dimensions of input (the dimension of the observation space) and output (the dimension of the principal subspace).
MultivariateStats.projection
— Methodprojection(M::AbstractMultispati)
Returns the projection matrix (of size (d, p)
). Each column of the projection matrix corresponds to a eigenvector. The eigenvectors are arranged in ascending order of the eigenvalues.
LinearAlgebra.eigvecs
— Methodeigvecs(M::AbstractMultispati)
Get the eigenvectors of the model M
.
LinearAlgebra.eigvals
— Methodeigvals(M::AbstractMultispati)
Get the eigenvalues of the model M
.
spatialPCA
Multispati.SpatialPCA
— TypeSpatialPCA <: AbstractMultispati
StatsAPI.fit
— Methodfit(SpatialPCA, X, W; ...)
Perform SpatialPCA over the data given a matrix X
and W
. Each column of X
is an observation. W
is a connectivity matrix where $w_{ij}$ is the connection from j -> i.
Keyword arguments
maxoutdim
: The output dimension, i.e. dimension of the transformed space (min(d, nc-1))solver
: The choice of solver::eig
: usesLinearAlgebra.eigen
(default):eigs
: usesArpack.eigs
(always used for sparse data)
tol
: Convergence tolerance foreigs
solver (default0.0
)maxiter
: Maximum number of iterations foreigs
solver (default300
)center_sparse
: Center sparse matrixX
(dense X will always be centered) (defaultfalse
)
References
StatsAPI.predict
— Methodpredict(M::SpatialPCA, x::AbstractVecOrMat{<:Real})
Transform the observations x
with the SpatialPCA model M
.
Here, x
can be either a vector of length d
or a matrix where each column is an observation.
MultivariateStats.reconstruct
— Methodreconstruct(M::SpatialPCA, y::AbstractVecOrMat{<:Real})
Approximately reconstruct the observations y
to the original space using the SpatialPCA model M
.
Here, y
can be either a vector of length p
or a matrix where each column gives the components for an observation.
Multispati.varianceMoransIdecomposition
— FunctionvarianceMoransIdecomposition(M::SpatialPCA, X)
Decompose the eigenvalues into a variance and Moran's I contribution given the model M
and matrix X
which was used for fitting the model.
Statistics.mean
— Methodmean(M::SpatialPCA)
Returns the mean vector (of length d
).