| AL_RNAmf | R Documentation |
The function acquires the new point and fidelity level by maximizing one of the four active learning criteria: ALM, ALC, ALD, or ALMC.
ALM (Active Learning MacKay): 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.
ALD (Active Learning Decomposition): It calculates the ALD criterion \frac{V_l(\bm{x})}{\sum^l_{j=1}C_j},
where V_l(\bm{x}) is the variance contribution of GP emulator
at each fidelity level l and C_j is the simulation cost at level j.
ALC (Active Learning Cohn): It calculates the ALC criterion
\frac{\Delta \sigma_L^{2}(l,\bm{x})}{\sum^l_{j=1}C_j} =
\frac{\int_{\Omega} \sigma_L^{*2}(\bm{\xi})-\tilde{\sigma}_L^{*2}(\bm{\xi};l,\bm{x}){\rm{d}}\bm{\xi}}{\sum^l_{j=1}C_j},
where f_L is the highest-fidelity simulation code,
\sigma_L^{*2}(\bm{\xi}) is the posterior variance of f_L(\bm{\xi}),
C_j is the simulation cost at fidelity level j,
and \tilde{\sigma}_L^{*2}(\bm{\xi};l,\bm{x}) is the posterior variance
based on the augmented design combining the current design and a new input location \bm{x}
at each fidelity level lower than or equal to l.
The integration is approximated by MC integration using uniform reference samples.
ALMC (Active Learning MacKay-Cohn): A hybrid approach.
It finds the optimal input location \bm{x}^*
by maximizing \sigma^{*2}_L(\bm{x}),
the posterior predictive variance at the highest-fidelity level L.
After selecting \bm{x}^*,
it finds the optimal fidelity level by maximizing ALC criterion at \bm{x}^*,
\text{argmax}_{l\in\{1,\ldots,L\}} \frac{\Delta \sigma_L^{2}(l,\bm{x}^*)}{\sum^l_{j=1}C_j},
where C_j is the simulation cost at level j.
A new point is acquired on Xcand. If Xcand=NULL and Xref=NULL, a new point is acquired on unit hypercube [0,1]^d.
For details, see Heo and Sung (2025, <\Sexpr[results=rd]{tools:::Rd_expr_doi("https://doi.org/10.1080/00401706.2024.2376173")}>).
AL_RNAmf(criterion = c("ALM", "ALC", "ALD", "ALMC"), fit,
Xref = NULL, Xcand = NULL, MC = FALSE, mc.sample = 100, cost = NULL,
use_optim = TRUE, parallel = FALSE, ncore = 1, trace = TRUE)
criterion |
character string specifying the active learning criterion to use. Must be one of |
fit |
object of class |
Xref |
vector or matrix of reference locations to approximate the integral of ALC. If |
Xcand |
vector or matrix of the candidate set for grid-based search. If |
MC |
logical indicating whether to use Monte Carlo approximation to impute the posterior variance (for ALC/ALMC). If |
mc.sample |
a number of MC samples generated for the imputation through MC approximation. Default is |
cost |
vector of the costs for each level of fidelity. If |
use_optim |
logical indicating whether to optimize the criterion using |
parallel |
logical indicating whether to use parallel computation. Default is |
ncore |
integer specifying the number of cores for parallel computation. Used only if |
trace |
logical indicating whether to print the computational progress and time. Default is |
For "ALC", or "ALMC", Xref plays a role of \bm{\xi} to approximate the integration.
To impute the posterior variance based on the augmented design \tilde{\sigma}_L^{*2}(\bm{\xi};l,\bm{x}),
MC approximation is used.
Due to the nested assumption, imputing y^{[s]}_{n_s+1} for each 1\leq s\leq l by drawing samples
from the posterior distribution of f_s(\bm{x}^{[s]}_{n_s+1})
based on the current design allows to compute \tilde{\sigma}_L^{*2}(\bm{\xi};l,\bm{x}).
Inverse of covariance matrix is computed by the Sherman-Morrison formula.
To search for the next acquisition \bm{x}^* by maximizing AL criterion,
the gradient-based optimization can be used by optim=TRUE.
Firstly, \tilde{\sigma}_L^{*2}(\bm{\xi};l,\bm{x}) is computed at
5 \times d number of points.
After that, the point minimizing \tilde{\sigma}_L^{*2}(\bm{\xi};l,\bm{x})
serves as a starting point of optimization by L-BFGS-B method.
Otherwise, when optim=FALSE, AL criterion is optimized only on Xcand.
The point is selected by maximizing the ALC criterion:
\text{argmax}_{l\in\{1,\ldots,L\}; \bm{x} \in \Omega}
\frac{\Delta \sigma_L^{2}(l,\bm{x})}{\sum^l_{j=1}C_j}.
A list containing:
AL: The values of the selected criterion at the candidate points (if use_optim=FALSE) or the optimized value (if use_optim=TRUE). For ALMC, it returns the ALC scores for each level at the chosen point.
cost: A copy of the cost argument.
Xcand: A copy of the Xcand argument used.
chosen: A list containing the chosen fidelity level and the new input location Xnext.
time: The computation time in seconds.
### 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]]
y1 <- f1(X1)
y2 <- f2(X2)
### n=100 uniform test data ###
x <- seq(0, 1, length.out = 100)
### fit an RNAmf ###
fit.RNAmf <- RNAmf(list(X1, X2), list(y1, y2), kernel = "sqex", constant=TRUE)
### 1. ALM Criterion ###
alm.RNAmf <- AL_RNAmf(criterion="ALM",
Xcand = x, fit=fit.RNAmf, cost = cost,
use_optim = FALSE, parallel = TRUE, ncore = 2)
print(alm.RNAmf$chosen)
### visualize ALM ###
oldpar <- par(mfrow = c(1, 2))
plot(x, alm.RNAmf$AL$ALM1,
type = "l", lty = 2,
xlab = "x", ylab = "ALM criterion at the low-fidelity level",
ylim = c(min(c(alm.RNAmf$AL$ALM1, alm.RNAmf$AL$ALM2)),
max(c(alm.RNAmf$AL$ALM1, alm.RNAmf$AL$ALM2))))
points(alm.RNAmf$chosen$Xnext,
alm.RNAmf$AL$ALM1[which(x == drop(alm.RNAmf$chosen$Xnext))],
pch = 16, cex = 1, col = "red")
plot(x, alm.RNAmf$AL$ALM2,
type = "l", lty = 2,
xlab = "x", ylab = "ALM criterion at the high-fidelity level",
ylim = c(min(c(alm.RNAmf$AL$ALM1, alm.RNAmf$AL$ALM2)),
max(c(alm.RNAmf$AL$ALM1, alm.RNAmf$AL$ALM2))))
par(oldpar)
### 2. ALD Criterion ###
ald.RNAmf <- AL_RNAmf(criterion="ALD",
Xcand = x, fit=fit.RNAmf, cost = cost,
use_optim = FALSE, parallel = TRUE, ncore = 2)
print(ald.RNAmf$chosen)
### visualize ALD ###
oldpar <- par(mfrow = c(1, 2))
plot(x, ald.RNAmf$AL$ALD1,
type = "l", lty = 2,
xlab = "x", ylab = "ALD criterion at the low-fidelity level",
ylim = c(min(c(ald.RNAmf$AL$ALD1, ald.RNAmf$AL$ALD2)),
max(c(ald.RNAmf$AL$ALD1, ald.RNAmf$AL$ALD2))))
points(ald.RNAmf$chosen$Xnext,
ald.RNAmf$AL$ALD1[which(x == drop(ald.RNAmf$chosen$Xnext))],
pch = 16, cex = 1, col = "red")
plot(x, ald.RNAmf$AL$ALD2,
type = "l", lty = 2,
xlab = "x", ylab = "ALD criterion at the high-fidelity level",
ylim = c(min(c(ald.RNAmf$AL$ALD1, ald.RNAmf$AL$ALD2)),
max(c(ald.RNAmf$AL$ALD1, ald.RNAmf$AL$ALD2))))
par(oldpar)
### 3. ALC Criterion ###
alc.RNAmf <- AL_RNAmf(criterion="ALC",
Xref = x, Xcand = x, fit=fit.RNAmf, cost = cost,
use_optim = FALSE, parallel = TRUE, ncore = 2)
print(alc.RNAmf$chosen)
### visualize ALC ###
oldpar <- par(mfrow = c(1, 2))
plot(x, alc.RNAmf$AL$ALC1,
type = "l", lty = 2,
xlab = "x", ylab = "ALC criterion augmented at the low-fidelity level",
ylim = c(min(c(alc.RNAmf$AL$ALC1, alc.RNAmf$AL$ALC2)),
max(c(alc.RNAmf$AL$ALC1, alc.RNAmf$AL$ALC2))))
points(alc.RNAmf$chosen$Xnext,
alc.RNAmf$AL$ALC1[which(x == drop(alc.RNAmf$chosen$Xnext))],
pch = 16, cex = 1, col = "red")
plot(x, alc.RNAmf$AL$ALC2,
type = "l", lty = 2,
xlab = "x", ylab = "ALC criterion augmented at the high-fidelity level",
ylim = c(min(c(alc.RNAmf$AL$ALC1, alc.RNAmf$AL$ALC2)),
max(c(alc.RNAmf$AL$ALC1, alc.RNAmf$AL$ALC2))))
par(oldpar)
### 4. ALMC Criterion ###
almc.RNAmf <- AL_RNAmf(criterion="ALMC",
Xref = x, Xcand = x, fit=fit.RNAmf, cost = cost,
use_optim = FALSE, parallel = TRUE, ncore = 2)
print(almc.RNAmf$chosen)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.