ALC_RNAmf: find the next point by ALC criterion

View source: R/ALC.R

ALC_RNAmfR Documentation

find the next point by ALC criterion

Description

The function acquires the new point by the Active learning Cohn (ALC) criterion. 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.

A new point is acquired on Xcand. If Xcand=NULL and Xref=NULL, a new point is acquired on unit hypercube [0,1]^d.

Usage

ALC_RNAmf(Xref = NULL, Xcand = NULL, fit, MC = FALSE, mc.sample = 100,
cost = NULL, optim = TRUE, parallel = FALSE, ncore = 1, trace=TRUE)

Arguments

Xref

vector or matrix of reference location to approximate the integral of ALC. If Xref=NULL, 100 \times d points from 0 to 1 are generated by Latin hypercube design. Default is NULL.

Xcand

vector or matrix of candidate set which could be added into the current design only when optim=FALSE. Xcand is the set of the points where ALC criterion is evaluated. If Xcand=NULL, Xref is used. Default is NULL. See details.

fit

object of class RNAmf.

MC

logical indicating whether to use MC approximation to impute the posterior variance or use posterior means. Default FALSE.

mc.sample

a number of mc samples generated for the imputation through MC approximation. Default is 100.

cost

vector of the costs for each level of fidelity. If cost=NULL, total costs at all levels would be 1. cost is encouraged to have a ascending order of positive value. Default is NULL.

optim

logical indicating whether to optimize AL criterion by optim's gradient-based L-BFGS-B method. If optim=TRUE, 5 \times d starting points are generated by Latin hypercube design for optimization. If optim=FALSE, AL criterion is optimized on the Xcand. Default is TRUE.

parallel

logical indicating whether to compute the AL criterion in parallel or not. If parallel=TRUE, parallel computation is utilized. Default is FALSE.

ncore

a number of core for parallel. It is only used if parallel=TRUE. Default is 1.

trace

logical indicating whether to print the computational time for each step. If trace=TRUE, the computation time for each step is printed. Default is TRUE.

Details

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. For details, see Heo and Sung (2025, <\Sexpr[results=rd]{tools:::Rd_expr_doi("https://doi.org/10.1080/00401706.2024.2376173")}>).

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 on the 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}.

Value

  • ALC: list of ALC criterion integrated on Xref when each data point on Xcand is added at each level l if optim=FALSE. If optim=TRUE, ALC returns NULL.

  • cost: a copy of cost.

  • Xcand: a copy of Xcand.

  • 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(list(X1, X2), list(y1, y2), kernel = "sqex", constant=TRUE)

### predict ###
predy <- predict(fit.RNAmf, x)$mu
predsig2 <- predict(fit.RNAmf, x)$sig2

### active learning with optim=TRUE ###
alc.RNAmf.optim <- ALC_RNAmf(
  Xref = x, Xcand = x, fit.RNAmf, cost = cost,
  optim = TRUE, parallel = TRUE, ncore = 2
)
print(alc.RNAmf.optim$time) # computation time of optim=TRUE

### active learning with optim=FALSE ###
alc.RNAmf <- ALC_RNAmf(
  Xref = x, Xcand = x, fit.RNAmf, cost = cost,
  optim = FALSE, parallel = TRUE, ncore = 2
)
print(alc.RNAmf$time) # computation time of optim=FALSE

### visualize ALC ###
oldpar <- par(mfrow = c(1, 2))
plot(x, alc.RNAmf$ALC$ALC1,
  type = "l", lty = 2,
  xlab = "x", ylab = "ALC criterion augmented at the low-fidelity level",
  ylim = c(min(c(alc.RNAmf$ALC$ALC1, alc.RNAmf$ALC$ALC2)),
           max(c(alc.RNAmf$ALC$ALC1, alc.RNAmf$ALC$ALC2)))
)
points(alc.RNAmf$chosen$Xnext,
  alc.RNAmf$ALC$ALC1[which(x == drop(alc.RNAmf$chosen$Xnext))],
  pch = 16, cex = 1, col = "red"
)
plot(x, alc.RNAmf$ALC$ALC2,
  type = "l", lty = 2,
  xlab = "x", ylab = "ALC criterion augmented at the high-fidelity level",
  ylim = c(min(c(alc.RNAmf$ALC$ALC1, alc.RNAmf$ALC$ALC2)),
           max(c(alc.RNAmf$ALC$ALC1, alc.RNAmf$ALC$ALC2)))
)
par(oldpar)


RNAmf documentation built on Sept. 13, 2025, 1:09 a.m.

Related to ALC_RNAmf in RNAmf...