ALMC_RNAmf: find the next point by ALMC criterion

View source: R/ALMC.R

ALMC_RNAmfR Documentation

find the next point by ALMC criterion

Description

The function acquires the new point by the hybrid approach, referred to as Active learning MacKay-Cohn (ALMC) criterion. 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. See ALC_RNAmf. For details, see Heo and Sung (2025, <\Sexpr[results=rd]{tools:::Rd_expr_doi("https://doi.org/10.1080/00401706.2024.2376173")}>).

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

ALMC_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 ALM criterion is evaluated. If Xcand=NULL, Xref is used. Default is NULL.

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.

Value

  • ALMC: vector of ALMC criterion \frac{\Delta \sigma_L^{2}(l,\bm{x^*})}{\sum^l_{j=1}C_j} for 1\leq l\leq L.

  • ALM: vector of ALM criterion computed at each point of Xcand at the highest fidelity level if optim=FALSE. If optim=TRUE, ALM returns NULL.

  • 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 ###
almc.RNAmf.optim <- ALMC_RNAmf(
  Xref = x, Xcand = x, fit.RNAmf, cost = cost,
  optim = TRUE, parallel = TRUE, ncore = 2
)
print(almc.RNAmf.optim$time) # computation time of optim=TRUE

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

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


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

Related to ALMC_RNAmf in RNAmf...