IMSPE_AL | R Documentation |
Search for the best next design point according to the IMSPE
criterion given a current MuFiMeshGP
model fit.
IMSPE_AL(
object,
t.min,
t.max,
cost.func,
cost.new = 0,
gr = FALSE,
gr_cost.func = NULL,
DesCand = NULL,
Wijs = NULL,
Hijs = NULL,
control = list(multi.start.n = 20, maxit = 20, DesStart = NULL, seed = NULL, ncores =
1)
)
object |
Current |
t.min , t.max |
Lower and upper bounds on the fidelity space for the search. |
cost.func |
Function that maps the tunable parameter |
cost.new |
(optional) Cost of running a new simulation at a new fidelity level, scalar. |
gr |
whether the gradient should be used in the optimization of the IMSPE. (Not recommended due to numerical errors) |
gr_cost.func |
If |
DesCand |
Design candidates to evaluate from. |
Wijs , Hijs |
(optional) Matrices from previous IMSPE search to obtain faster computation through matrix decomposition. |
control |
list of arguments udes for the optimization. |
a list with:
x
: the optimal input parameter, a vector.
t
: the optimal tuning parameter, a scalar.
value
: the IMSPE reduction at the optimal design location.
new
: whether the optimal tuning parameter defines a new fidelity level.
id
: the index of the optimal design location if DesCand
is used.
# Example code
f <- function(x, t){
x <- c(x)
return(exp(-1.4*x)*cos(3.5*pi*x)+sin(40*x)/10*t^2)
}
set.seed(1)
X <- matrix(runif(15,0,1), ncol = 1)
tt <- runif(15,0.5,2)
Y <- f(c(X), tt)
fit.mufimeshgp <- MuFiMeshGP(X, tt, Y)
xx <- matrix(seq(0,1,0.01), ncol = 1)
ftrue <- f(xx, 0)
# predict
pred.mufimeshgp <- predict(fit.mufimeshgp, xx, rep(0,101))
mu <- pred.mufimeshgp$mean
s <- pred.mufimeshgp$sd
lower <- mu + qnorm(0.025)*s
upper <- mu + qnorm(0.975)*s
# plot
oldpar <- par(mfrow = c(1,2))
plot(xx, ftrue, "l", ylim = c(-1,1.3), ylab = "y", xlab = "x")
lines(c(xx), mu, col = "blue")
lines(c(xx), lower, col = "blue", lty = 2)
lines(c(xx), upper, col = "blue", lty = 2)
points(c(X), Y, col = "red")
### RMSE ###
print(sqrt(mean((ftrue - mu))^2))
best <- IMSPE_AL(fit.mufimeshgp, 0.5, 2, function(t) return(1 / t^2))
new.Y <- f(best$x, best$t)
fit.mufimeshgp <- update(fit.mufimeshgp, best$x, best$t, new.Y)
pred.mufimeshgp <- predict(fit.mufimeshgp, xx, rep(0, 101))
mu <- pred.mufimeshgp$mean
s <- pred.mufimeshgp$sd
lower <- mu + qnorm(0.025)*s
upper <- mu + qnorm(0.975)*s
plot(xx, ftrue, "l", ylim = c(-1,1.3), ylab = "y", xlab = "x")
lines(c(xx), mu, col = "blue")
lines(c(xx), lower, col = "blue", lty = 2)
lines(c(xx), upper, col = "blue", lty = 2)
points(c(X), Y, col = "red")
points(c(best$x), new.Y, col = "green")
par(oldpar)
### RMSE ###
print(sqrt(mean((ftrue - mu))^2))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.