LatentFactorKernel {GauPro} | R Documentation |
Latent Factor Kernel R6 class
Description
Latent Factor Kernel R6 class
Latent Factor Kernel R6 class
Usage
k_LatentFactorKernel(
s2 = 1,
D,
nlevels,
xindex,
latentdim,
p_lower = 0,
p_upper = 1,
p_est = TRUE,
s2_lower = 1e-08,
s2_upper = 1e+08,
s2_est = TRUE,
useC = TRUE,
offdiagequal = 1 - 1e-06
)
Arguments
s2 |
Initial variance |
D |
Number of input dimensions of data |
nlevels |
Number of levels for the factor |
xindex |
Index of X to use the kernel on |
latentdim |
Dimension of embedding space |
p_lower |
Lower bound for p |
p_upper |
Upper bound for p |
p_est |
Should p be estimated? |
s2_lower |
Lower bound for s2 |
s2_upper |
Upper bound for s2 |
s2_est |
Should s2 be estimated? |
useC |
Should C code used? Much faster. |
offdiagequal |
What should offdiagonal values be set to when the indices are the same? Use to avoid decomposition errors, similar to adding a nugget. |
Format
R6Class
object.
Details
Used for factor variables, a single dimension. Each level of the factor gets mapped into a latent space, then the distances in that space determine their correlations.
Value
Object of R6Class
with methods for fitting GP model.
Super class
GauPro::GauPro_kernel
-> GauPro_kernel_LatentFactorKernel
Public fields
p
Parameter for correlation
p_est
Should p be estimated?
p_lower
Lower bound of p
p_upper
Upper bound of p
p_length
length of p
s2
variance
s2_est
Is s2 estimated?
logs2
Log of s2
logs2_lower
Lower bound of logs2
logs2_upper
Upper bound of logs2
xindex
Index of the factor (which column of X)
nlevels
Number of levels for the factor
latentdim
Dimension of embedding space
pf_to_p_log
Logical vector used to convert pf to p
p_to_pf_inds
Vector of indexes used to convert p to pf
offdiagequal
What should offdiagonal values be set to when the indices are the same? Use to avoid decomposition errors, similar to adding a nugget.
Methods
Public methods
Inherited methods
Method new()
Initialize kernel object
Usage
LatentFactorKernel$new( s2 = 1, D, nlevels, xindex, latentdim, p_lower = 0, p_upper = 1, p_est = TRUE, s2_lower = 1e-08, s2_upper = 1e+08, s2_est = TRUE, useC = TRUE, offdiagequal = 1 - 1e-06 )
Arguments
s2
Initial variance
D
Number of input dimensions of data
nlevels
Number of levels for the factor
xindex
Index of X to use the kernel on
latentdim
Dimension of embedding space
p_lower
Lower bound for p
p_upper
Upper bound for p
p_est
Should p be estimated?
s2_lower
Lower bound for s2
s2_upper
Upper bound for s2
s2_est
Should s2 be estimated?
useC
Should C code used? Much faster.
offdiagequal
What should offdiagonal values be set to when the indices are the same? Use to avoid decomposition errors, similar to adding a nugget.
Method k()
Calculate covariance between two points
Usage
LatentFactorKernel$k(x, y = NULL, p = self$p, s2 = self$s2, params = NULL)
Arguments
x
vector.
y
vector, optional. If excluded, find correlation of x with itself.
p
Correlation parameters.
s2
Variance parameter.
params
parameters to use instead of beta and s2.
Method kone()
Find covariance of two points
Usage
LatentFactorKernel$kone( x, y, pf, s2, isdiag = TRUE, offdiagequal = self$offdiagequal )
Arguments
x
vector
y
vector
pf
correlation parameters on regular scale, includes zeroes for first level.
s2
Variance parameter
isdiag
Is this on the diagonal of the covariance?
offdiagequal
What should offdiagonal values be set to when the indices are the same? Use to avoid decomposition errors, similar to adding a nugget.
Method dC_dparams()
Derivative of covariance with respect to parameters
Usage
LatentFactorKernel$dC_dparams(params = NULL, X, C_nonug, C, nug)
Arguments
params
Kernel parameters
X
matrix of points in rows
C_nonug
Covariance without nugget added to diagonal
C
Covariance with nugget
nug
Value of nugget
Method C_dC_dparams()
Calculate covariance matrix and its derivative with respect to parameters
Usage
LatentFactorKernel$C_dC_dparams(params = NULL, X, nug)
Arguments
params
Kernel parameters
X
matrix of points in rows
nug
Value of nugget
Method dC_dx()
Derivative of covariance with respect to X
Usage
LatentFactorKernel$dC_dx(XX, X, ...)
Arguments
XX
matrix of points
X
matrix of points to take derivative with respect to
...
Additional args, not used
Method param_optim_start()
Starting point for parameters for optimization
Usage
LatentFactorKernel$param_optim_start( jitter = F, y, p_est = self$p_est, s2_est = self$s2_est )
Arguments
jitter
Should there be a jitter?
y
Output
p_est
Is p being estimated?
s2_est
Is s2 being estimated?
Method param_optim_start0()
Starting point for parameters for optimization
Usage
LatentFactorKernel$param_optim_start0( jitter = F, y, p_est = self$p_est, s2_est = self$s2_est )
Arguments
jitter
Should there be a jitter?
y
Output
p_est
Is p being estimated?
s2_est
Is s2 being estimated?
Method param_optim_lower()
Lower bounds of parameters for optimization
Usage
LatentFactorKernel$param_optim_lower(p_est = self$p_est, s2_est = self$s2_est)
Arguments
p_est
Is p being estimated?
s2_est
Is s2 being estimated?
Method param_optim_upper()
Upper bounds of parameters for optimization
Usage
LatentFactorKernel$param_optim_upper(p_est = self$p_est, s2_est = self$s2_est)
Arguments
p_est
Is p being estimated?
s2_est
Is s2 being estimated?
Method set_params_from_optim()
Set parameters from optimization output
Usage
LatentFactorKernel$set_params_from_optim( optim_out, p_est = self$p_est, s2_est = self$s2_est )
Arguments
optim_out
Output from optimization
p_est
Is p being estimated?
s2_est
Is s2 being estimated?
Method p_to_pf()
Convert p (short parameter vector) to pf (long parameter vector with zeros).
Usage
LatentFactorKernel$p_to_pf(p)
Arguments
p
Parameter vector
Method s2_from_params()
Get s2 from params vector
Usage
LatentFactorKernel$s2_from_params(params, s2_est = self$s2_est)
Arguments
params
parameter vector
s2_est
Is s2 being estimated?
Method plotLatent()
Plot the points in the latent space
Usage
LatentFactorKernel$plotLatent()
Method print()
Print this object
Usage
LatentFactorKernel$print()
Method clone()
The objects of this class are cloneable with this method.
Usage
LatentFactorKernel$clone(deep = FALSE)
Arguments
deep
Whether to make a deep clone.
References
https://stackoverflow.com/questions/27086195/linear-index-upper-triangular-matrix
Examples
# Create a new kernel for a single factor with 5 levels,
# mapped into two latent dimensions.
kk <- LatentFactorKernel$new(D=1, nlevels=5, xindex=1, latentdim=2)
# Random initial parameter values
kk$p
# Plots to understand
kk$plotLatent()
kk$plot()
# 5 levels, 1/4 are similar and 2/3/5 are similar
n <- 30
x <- matrix(sample(1:5, n, TRUE))
y <- c(ifelse(x == 1 | x == 4, 4, -3) + rnorm(n,0,.1))
plot(c(x), y)
m5 <- GauPro_kernel_model$new(
X=x, Z=y,
kernel=LatentFactorKernel$new(D=1, nlevels = 5, xindex = 1, latentdim = 2))
m5$kernel$p
# We should see 1/4 and 2/3/4 in separate clusters
m5$kernel$plotLatent()
library(dplyr)
n <- 20
X <- cbind(matrix(runif(n,2,6), ncol=1),
matrix(sample(1:2, size=n, replace=TRUE), ncol=1))
X <- rbind(X, c(3.3,3), c(3.7,3))
n <- nrow(X)
Z <- X[,1] - (4-X[,2])^2 + rnorm(n,0,.1)
plot(X[,1], Z, col=X[,2])
tibble(X=X, Z) %>% arrange(X,Z)
k2a <- IgnoreIndsKernel$new(k=Gaussian$new(D=1), ignoreinds = 2)
k2b <- LatentFactorKernel$new(D=2, nlevels=3, xind=2, latentdim=2)
k2 <- k2a * k2b
k2b$p_upper <- .65*k2b$p_upper
gp <- GauPro_kernel_model$new(X=X, Z=Z, kernel = k2, verbose = 5,
nug.min=1e-2, restarts=1)
gp$kernel$k1$kernel$beta
gp$kernel$k2$p
gp$kernel$k(x = gp$X)
tibble(X=X, Z=Z, pred=gp$predict(X)) %>% arrange(X, Z)
tibble(X=X[,2], Z) %>% group_by(X) %>% summarize(n=n(), mean(Z))
curve(gp$pred(cbind(matrix(x,ncol=1),1)),2,6, ylim=c(min(Z), max(Z)))
points(X[X[,2]==1,1], Z[X[,2]==1])
curve(gp$pred(cbind(matrix(x,ncol=1),2)), add=TRUE, col=2)
points(X[X[,2]==2,1], Z[X[,2]==2], col=2)
curve(gp$pred(cbind(matrix(x,ncol=1),3)), add=TRUE, col=3)
points(X[X[,2]==3,1], Z[X[,2]==3], col=3)
legend(legend=1:3, fill=1:3, x="topleft")
# See which points affect (5.5, 3 themost)
data.frame(X, cov=gp$kernel$k(X, c(5.5,3))) %>% arrange(-cov)
plot(k2b)