ALM_RNAmf {RNAmf} | R Documentation |
find the next point by ALM criterion
Description
The function acquires the new point by the Active learning MacKay (ALM) criterion.
It calculates the ALM criterion \frac{\sigma^{*2}_l(\bm{x})}{\sum^l_{j=1}C_j}
,
where \sigma^{*2}_l(\bm{x})
is the posterior predictive variance
at each fidelity level l
and C_j
is the simulation cost at level j
.
For details, see Heo and Sung (2024, <doi:10.1080/00401706.2024.2376173>).
A new point is acquired on Xcand
. If Xcand=NULL
, a new point is acquired on unit hypercube [0,1]^d
.
Usage
ALM_RNAmf(Xcand = NULL, fit, cost = NULL, optim = TRUE,
parallel = FALSE, ncore = 1, trace=TRUE)
Arguments
Xcand |
vector or matrix of candidate set which could be added into the current design only used when |
fit |
object of class |
cost |
vector of the costs for each level of fidelity. If |
optim |
logical indicating whether to optimize AL criterion by |
parallel |
logical indicating whether to compute the AL criterion in parallel or not. If |
ncore |
a number of core for parallel. It is only used if |
trace |
logical indicating whether to print the computational time for each step. If |
Value
-
ALM
: list of ALM criterion computed at each point ofXcand
at each level ifoptim=FALSE
. Ifoptim=TRUE
,ALM
returnsNULL
. -
cost
: a copy ofcost
. -
Xcand
: a copy ofXcand
. -
chosen
: list of chosen level and point. -
time
: a scalar of the time for the computation.
Examples
library(lhs)
library(doParallel)
library(foreach)
### simulation costs ###
cost <- c(1, 3)
### 1-d Perdikaris function in Perdikaris, et al. (2017) ###
# low-fidelity function
f1 <- function(x) {
sin(8 * pi * x)
}
# high-fidelity function
f2 <- function(x) {
(x - sqrt(2)) * (sin(8 * pi * x))^2
}
### training data ###
n1 <- 13
n2 <- 8
### fix seed to reproduce the result ###
set.seed(1)
### generate initial nested design ###
X <- NestedX(c(n1, n2), 1)
X1 <- X[[1]]
X2 <- X[[2]]
### n1 and n2 might be changed from NestedX ###
### assign n1 and n2 again ###
n1 <- nrow(X1)
n2 <- nrow(X2)
y1 <- f1(X1)
y2 <- f2(X2)
### n=100 uniform test data ###
x <- seq(0, 1, length.out = 100)
### fit an RNAmf ###
fit.RNAmf <- RNAmf_two_level(X1, y1, X2, y2, kernel = "sqex")
### predict ###
predy <- predict(fit.RNAmf, x)$mu
predsig2 <- predict(fit.RNAmf, x)$sig2
### active learning with optim=TRUE ###
alm.RNAmf.optim <- ALM_RNAmf(
Xcand = x, fit.RNAmf, cost = cost,
optim = TRUE, parallel = TRUE, ncore = 2
)
print(alm.RNAmf.optim$time) # computation time of optim=TRUE
### active learning with optim=FALSE ###
alm.RNAmf <- ALM_RNAmf(
Xcand = x, fit.RNAmf, cost = cost,
optim = FALSE, parallel = TRUE, ncore = 2
)
print(alm.RNAmf$time) # computation time of optim=FALSE
### visualize ALM ###
oldpar <- par(mfrow = c(1, 2))
plot(x, alm.RNAmf$ALM$ALM1,
type = "l", lty = 2,
xlab = "x", ylab = "ALM criterion at the low-fidelity level",
ylim = c(min(c(alm.RNAmf$ALM$ALM1, alm.RNAmf$ALM$ALM2)),
max(c(alm.RNAmf$ALM$ALM1, alm.RNAmf$ALM$ALM2)))
)
points(alm.RNAmf$chosen$Xnext,
alm.RNAmf$ALM$ALM1[which(x == drop(alm.RNAmf$chosen$Xnext))],
pch = 16, cex = 1, col = "red"
)
plot(x, alm.RNAmf$ALM$ALM2,
type = "l", lty = 2,
xlab = "x", ylab = "ALM criterion at the high-fidelity level",
ylim = c(min(c(alm.RNAmf$ALM$ALM1, alm.RNAmf$ALM$ALM2)),
max(c(alm.RNAmf$ALM$ALM1, alm.RNAmf$ALM$ALM2)))
)
par(oldpar)