diff --git a/DESCRIPTION b/DESCRIPTION
index a10fe9d..4eb3ac5 100644
--- a/DESCRIPTION
+++ b/DESCRIPTION
@@ -1,6 +1,6 @@
Package: lhs
Title: Latin Hypercube Samples
-Version: 1.3.0
+Version: 1.4.0
Authors@R:
person(given = "Rob",
family = "Carnell",
@@ -8,7 +8,7 @@ Authors@R:
email = "bertcarnell@gmail.com",
comment = c(ORCID = "0009-0009-0465-7564")
)
-Description: Provides a number of methods for creating and augmenting Latin Hypercube Samples and Orthogonal Array Latin Hypercube Samples.
+Description: Provides a number of methods for creating and augmenting Latin Hypercube Samples and Orthogonal Array Latin Hypercube Samples, including a unified interface for generating designs with transformed margins.
License: GPL-3
Encoding: UTF-8
Depends: R (>= 4.1.0)
diff --git a/NAMESPACE b/NAMESPACE
index 03a7d76..1bbf2ad 100644
--- a/NAMESPACE
+++ b/NAMESPACE
@@ -15,6 +15,7 @@ export(create_oalhs)
export(geneticLHS)
export(get_library_versions)
export(improvedLHS)
+export(lhs_design)
export(maximinLHS)
export(oa_to_oalhs)
export(optAugmentLHS)
@@ -28,6 +29,7 @@ export(qfactor)
export(qinteger)
export(randomLHS)
export(runifint)
+export(transform_lhs)
import(Rcpp)
importFrom(stats,dist)
importFrom(stats,na.exclude)
diff --git a/NEWS b/NEWS
index e6c9f56..d9b9565 100644
--- a/NEWS
+++ b/NEWS
@@ -140,3 +140,11 @@ Changes in version 1.3.0 (2026-04-19)
- Added ability to interrupt from long running C++ code processes
- Updated R dependency to 4.0.0 and C++ to 17
+
+Changes in version 1.4.0 (2026-06-06)
+
+- Added lhs_design, a unified interface that generates a Latin hypercube sample
+ and transforms each margin to a distribution, factor, or integer range,
+ returning a data.frame
+- Added transform_lhs to apply the same margin transformations to an existing
+ Latin hypercube sample
diff --git a/R/lhs_design.R b/R/lhs_design.R
new file mode 100644
index 0000000..773d747
--- /dev/null
+++ b/R/lhs_design.R
@@ -0,0 +1,215 @@
+# Copyright 2026 Robert Carnell
+
+#' Create a Latin hypercube design with transformed margins
+#'
+#' A unified, high-level interface that generates a Latin hypercube sample and
+#' transforms each margin to a user supplied distribution, factor, or integer
+#' range in a single call, returning a ready to use \code{data.frame}.
+#'
+#' @details \code{lhs_design} is a convenience wrapper around the design
+#' generators ([randomLHS()], [improvedLHS()], [maximinLHS()], [optimumLHS()],
+#' [geneticLHS()], and [create_oalhs()]) and the quantile transformations
+#' ([qfactor()], [qinteger()], [qdirichlet()]). The number of columns of the
+#' underlying Latin hypercube sample is determined from the \code{variables}
+#' specification, the sample is generated using the requested \code{type}, and
+#' each margin is transformed and assembled into a named \code{data.frame}.
+#'
+#' [transform_lhs()] applies the same transformations to an existing Latin
+#' hypercube sample (for example one created by [augmentLHS()], [slicedLHS()],
+#' or any matrix with values on [0,1]) without generating a new design.
+#'
+#' Each element of \code{variables} describes one variable and may be one of:
+#' \describe{
+#' \item{a function}{a function of a single numeric vector of values on [0,1]
+#' that returns the transformed column, e.g. \code{function(p) qnorm(p, 2, 0.5)}.
+#' It consumes a single column of the Latin hypercube sample.}
+#' \item{a distribution list}{a list with a character element \code{dist}
+#' naming a quantile function (e.g. \code{"qnorm"}) and any further named
+#' arguments passed to that function, e.g.
+#' \code{list(dist = "qnorm", mean = 2, sd = 0.5)}. It consumes a single
+#' column.}
+#' \item{a factor list}{a list with \code{type = "factor"}, a \code{levels}
+#' vector, and an optional logical \code{ordered} (default \code{FALSE}). It
+#' consumes a single column and is transformed with [qfactor()].}
+#' \item{an integer list}{a list with \code{type = "integer"} and integer
+#' \code{min} and \code{max} values. It consumes a single column and is
+#' transformed with [qinteger()].}
+#' \item{a dirichlet list}{a list with \code{type = "dirichlet"} and a numeric
+#' \code{alpha} vector. It consumes \code{length(alpha)} columns, transforms
+#' them with [qdirichlet()], and produces \code{length(alpha)} output columns
+#' named with the variable name followed by an index.}
+#' }
+#'
+#' @param n the number of rows or samples
+#' @param variables a named list describing the marginal transformation of each
+#' variable. See Details for the accepted specifications. The names become the
+#' columns of the returned \code{data.frame}.
+#' @param type the type of Latin hypercube sample to generate. One of
+#' \code{"random"}, \code{"improved"}, \code{"maximin"}, \code{"optimum"},
+#' \code{"genetic"}, or \code{"oalhs"}.
+#' @param ... further arguments passed to the underlying design generator, for
+#' example \code{dup} for \code{"improved"} or \code{maxSweeps} for
+#' \code{"optimum"}.
+#' @return a \code{data.frame} with one column per output variable and values
+#' transformed to the requested margins
+#' @export
+#'
+#' @author Rob Carnell
+#'
+#' @keywords design
+#'
+#' @seealso [transform_lhs()] to transform an existing design. [randomLHS()],
+#' [improvedLHS()], [maximinLHS()], [optimumLHS()], [geneticLHS()], and
+#' [create_oalhs()] for the underlying generators. [qfactor()], [qinteger()],
+#' and [qdirichlet()] for the underlying transformations.
+#'
+#' @examples
+#' set.seed(1234)
+#' D <- lhs_design(20, variables = list(
+#' x = list(dist = "qnorm", mean = 2, sd = 0.5),
+#' grp = list(type = "factor", levels = c("a", "b", "c")),
+#' count = list(type = "integer", min = 5L, max = 17L),
+#' z = function(p) qexp(p, rate = 2),
+#' comp = list(type = "dirichlet", alpha = c(2, 3, 4))
+#' ), type = "maximin")
+#' head(D)
+lhs_design <- function(n, variables,
+ type = c("random", "improved", "maximin", "optimum",
+ "genetic", "oalhs"),
+ ...)
+{
+ type <- match.arg(type)
+ if (length(n) != 1 || is.na(n) || is.infinite(n) || n != floor(n) || n < 1)
+ stop("n must be a positive integer")
+ .lhs_check_variables(variables)
+
+ k <- sum(vapply(variables, .lhs_spec_ncol, integer(1)))
+
+ lhs <- switch(type,
+ random = randomLHS(n, k, ...),
+ improved = improvedLHS(n, k, ...),
+ maximin = maximinLHS(n, k, ...),
+ optimum = optimumLHS(n, k, ...),
+ genetic = geneticLHS(n, k, ...),
+ oalhs = .lhs_call_oalhs(n, k, ...))
+
+ transform_lhs(lhs, variables)
+}
+
+#' @rdname lhs_design
+#'
+#' @param lhs a Latin hypercube sample matrix with values on [0,1]. The number
+#' of columns must equal the total number of columns consumed by
+#' \code{variables}.
+#' @export
+transform_lhs <- function(lhs, variables)
+{
+ if (!is.matrix(lhs))
+ stop("lhs must be a matrix")
+ if (any(is.na(lhs)))
+ stop("lhs cannot contain any NA entries")
+ if (any(lhs < 0 | lhs > 1))
+ stop("lhs must have entries on the interval [0,1]")
+ .lhs_check_variables(variables)
+
+ ncols <- vapply(variables, .lhs_spec_ncol, integer(1))
+ if (sum(ncols) != ncol(lhs))
+ stop(paste0("the number of columns of lhs (", ncol(lhs),
+ ") must equal the number of columns required by variables (",
+ sum(ncols), ")"))
+
+ out <- list()
+ col <- 1L
+ for (i in seq_along(variables))
+ {
+ nm <- names(variables)[i]
+ m <- ncols[i]
+ idx <- seq.int(col, length.out = m)
+ if (m == 1L)
+ {
+ out[[nm]] <- .lhs_apply_spec(lhs[, idx], variables[[i]])
+ } else
+ {
+ res <- .lhs_apply_spec(lhs[, idx, drop = FALSE], variables[[i]])
+ for (jj in seq_len(m))
+ out[[paste0(nm, jj)]] <- res[, jj]
+ }
+ col <- col + m
+ }
+
+ data.frame(out, stringsAsFactors = FALSE, check.names = FALSE)
+}
+
+# Validate the variables list shared by lhs_design and transform_lhs
+.lhs_check_variables <- function(variables)
+{
+ if (!is.list(variables) || length(variables) == 0)
+ stop("variables must be a non-empty list")
+ nm <- names(variables)
+ if (is.null(nm) || any(!nzchar(nm)))
+ stop("every element of variables must be named")
+ if (anyDuplicated(nm))
+ stop("the names of variables must be unique")
+ invisible(TRUE)
+}
+
+# Number of Latin hypercube columns consumed by a single variable specification
+.lhs_spec_ncol <- function(spec)
+{
+ if (is.function(spec))
+ return(1L)
+ if (is.list(spec))
+ {
+ if (identical(spec$type, "dirichlet"))
+ {
+ if (is.null(spec$alpha))
+ stop("a dirichlet specification must include 'alpha'")
+ return(length(spec$alpha))
+ }
+ return(1L)
+ }
+ stop("each variable specification must be a function or a list")
+}
+
+# Apply a single variable specification to one or more columns of an LHS
+.lhs_apply_spec <- function(x, spec)
+{
+ if (is.function(spec))
+ return(spec(x))
+
+ type <- spec$type
+ if (!is.null(type) && type == "factor")
+ {
+ if (is.null(spec$levels))
+ stop("a factor specification must include 'levels'")
+ fact <- factor(spec$levels, levels = spec$levels,
+ ordered = isTRUE(spec$ordered))
+ return(qfactor(x, fact))
+ }
+ if (!is.null(type) && type == "integer")
+ {
+ if (is.null(spec$min) || is.null(spec$max))
+ stop("an integer specification must include 'min' and 'max'")
+ return(qinteger(x, spec$min, spec$max))
+ }
+ if (!is.null(type) && type == "dirichlet")
+ return(qdirichlet(x, spec$alpha))
+ if (!is.null(spec$dist))
+ {
+ fn <- match.fun(spec$dist)
+ args <- spec[setdiff(names(spec), "dist")]
+ return(do.call(fn, c(list(x), args)))
+ }
+ stop("unrecognized variable specification")
+}
+
+# Call create_oalhs supplying defaults for its required flags
+.lhs_call_oalhs <- function(n, k, ...)
+{
+ args <- list(...)
+ if (is.null(args$bChooseLargerDesign))
+ args$bChooseLargerDesign <- TRUE
+ if (is.null(args$bverbose))
+ args$bverbose <- FALSE
+ do.call(create_oalhs, c(list(n, k), args))
+}
diff --git a/README.Rmd b/README.Rmd
index 037bed2..f77c9a2 100644
--- a/README.Rmd
+++ b/README.Rmd
@@ -131,6 +131,20 @@ W16 <- create_oalhs(10, 3, bChooseLargerDesign = TRUE, bverbose = FALSE)
dim(W16)
```
+Generate a design and transform each margin to a distribution, factor, or
+integer range in a single call:
+
+```{r design, echo=TRUE}
+set.seed(42)
+D <- lhs_design(8, variables = list(
+ temperature = list(dist = "qnorm", mean = 350, sd = 10),
+ material = list(type = "factor", levels = c("steel", "alloy", "ti")),
+ cycles = list(type = "integer", min = 1L, max = 100L)
+), type = "maximin")
+dim(D)
+sapply(D, class)
+```
+
## Help
R-Help Examples of using the LHS package
diff --git a/README.md b/README.md
index c130de1..d0d3554 100644
--- a/README.md
+++ b/README.md
@@ -1,21 +1,26 @@
+
+
+
|
+
+
-| Actions | Code Coverage | Website | Doxygen | CRAN Downloads | CRAN |
-|:-----------------------------------------------------------------------------------------------------------------------------------:|:--------------------------------------------------------------------------------------------------------------------------------------------------:|:--------------------------------------------------------------------------------------------:|:-----------------------------------------------------------------------------------------------------------:|:------------------------------------------------------------------------------------:|:--------------------------------------------------------------------------------------------------:|
+| Actions | Code Coverage | Website | Doxygen | CRAN Downloads | CRAN |
+|:--:|:--:|:--:|:--:|:--:|:--:|
| [](https://github.com/bertcarnell/lhs/actions) | [](https://codecov.io/github/bertcarnell/lhs?branch=master) | [](https://bertcarnell.github.io/lhs/) | [](https://bertcarnell.github.io/lhs/html/index.html) | [](https://cran.r-project.org/package=lhs) | [](https://cran.r-project.org/package=lhs) |
-| R-Universe | CRAN Packages | DOI |
-|:------------------------------------------------------------------------------------------------------------:|:---------------------------------------------------------------------------------------------------------------------------------------------------------:|:----------------------------------------------------------------------------------------------------------------------------:|
+| R-Universe | CRAN Packages | DOI |
+|:--:|:--:|:--:|
| [](https://bertcarnell.r-universe.dev/) | [](https://cran.r-project.org/web/checks/check_results_bertcarnell_at_gmail.com.html) | [](https://cran.r-project.org/package=lhs) |
# lhs
@@ -110,6 +115,28 @@ dim(W16)
## [1] 16 3
+Generate a design and transform each margin to a distribution, factor,
+or integer range in a single call:
+
+``` r
+set.seed(42)
+D <- lhs_design(8, variables = list(
+ temperature = list(dist = "qnorm", mean = 350, sd = 10),
+ material = list(type = "factor", levels = c("steel", "alloy", "ti")),
+ cycles = list(type = "integer", min = 1L, max = 100L)
+), type = "maximin")
+dim(D)
+```
+
+ ## [1] 8 3
+
+``` r
+sapply(D, class)
+```
+
+ ## temperature material cycles
+ ## "numeric" "factor" "numeric"
+
## Help
R-Help Examples of using the LHS package
diff --git a/man/lhs_design.Rd b/man/lhs_design.Rd
new file mode 100644
index 0000000..0f67766
--- /dev/null
+++ b/man/lhs_design.Rd
@@ -0,0 +1,100 @@
+% Generated by roxygen2: do not edit by hand
+% Please edit documentation in R/lhs_design.R
+\name{lhs_design}
+\alias{lhs_design}
+\alias{transform_lhs}
+\title{Create a Latin hypercube design with transformed margins}
+\usage{
+lhs_design(
+ n,
+ variables,
+ type = c("random", "improved", "maximin", "optimum", "genetic", "oalhs"),
+ ...
+)
+
+transform_lhs(lhs, variables)
+}
+\arguments{
+\item{n}{the number of rows or samples}
+
+\item{variables}{a named list describing the marginal transformation of each
+variable. See Details for the accepted specifications. The names become the
+columns of the returned \code{data.frame}.}
+
+\item{type}{the type of Latin hypercube sample to generate. One of
+\code{"random"}, \code{"improved"}, \code{"maximin"}, \code{"optimum"},
+\code{"genetic"}, or \code{"oalhs"}.}
+
+\item{...}{further arguments passed to the underlying design generator, for
+example \code{dup} for \code{"improved"} or \code{maxSweeps} for
+\code{"optimum"}.}
+
+\item{lhs}{a Latin hypercube sample matrix with values on [0,1]. The number
+of columns must equal the total number of columns consumed by
+\code{variables}.}
+}
+\value{
+a \code{data.frame} with one column per output variable and values
+transformed to the requested margins
+}
+\description{
+A unified, high-level interface that generates a Latin hypercube sample and
+transforms each margin to a user supplied distribution, factor, or integer
+range in a single call, returning a ready to use \code{data.frame}.
+}
+\details{
+\code{lhs_design} is a convenience wrapper around the design
+generators ([randomLHS()], [improvedLHS()], [maximinLHS()], [optimumLHS()],
+[geneticLHS()], and [create_oalhs()]) and the quantile transformations
+([qfactor()], [qinteger()], [qdirichlet()]). The number of columns of the
+underlying Latin hypercube sample is determined from the \code{variables}
+specification, the sample is generated using the requested \code{type}, and
+each margin is transformed and assembled into a named \code{data.frame}.
+
+[transform_lhs()] applies the same transformations to an existing Latin
+hypercube sample (for example one created by [augmentLHS()], [slicedLHS()],
+or any matrix with values on [0,1]) without generating a new design.
+
+Each element of \code{variables} describes one variable and may be one of:
+\describe{
+ \item{a function}{a function of a single numeric vector of values on [0,1]
+ that returns the transformed column, e.g. \code{function(p) qnorm(p, 2, 0.5)}.
+ It consumes a single column of the Latin hypercube sample.}
+ \item{a distribution list}{a list with a character element \code{dist}
+ naming a quantile function (e.g. \code{"qnorm"}) and any further named
+ arguments passed to that function, e.g.
+ \code{list(dist = "qnorm", mean = 2, sd = 0.5)}. It consumes a single
+ column.}
+ \item{a factor list}{a list with \code{type = "factor"}, a \code{levels}
+ vector, and an optional logical \code{ordered} (default \code{FALSE}). It
+ consumes a single column and is transformed with [qfactor()].}
+ \item{an integer list}{a list with \code{type = "integer"} and integer
+ \code{min} and \code{max} values. It consumes a single column and is
+ transformed with [qinteger()].}
+ \item{a dirichlet list}{a list with \code{type = "dirichlet"} and a numeric
+ \code{alpha} vector. It consumes \code{length(alpha)} columns, transforms
+ them with [qdirichlet()], and produces \code{length(alpha)} output columns
+ named with the variable name followed by an index.}
+}
+}
+\examples{
+set.seed(1234)
+D <- lhs_design(20, variables = list(
+ x = list(dist = "qnorm", mean = 2, sd = 0.5),
+ grp = list(type = "factor", levels = c("a", "b", "c")),
+ count = list(type = "integer", min = 5L, max = 17L),
+ z = function(p) qexp(p, rate = 2),
+ comp = list(type = "dirichlet", alpha = c(2, 3, 4))
+), type = "maximin")
+head(D)
+}
+\seealso{
+[transform_lhs()] to transform an existing design. [randomLHS()],
+[improvedLHS()], [maximinLHS()], [optimumLHS()], [geneticLHS()], and
+[create_oalhs()] for the underlying generators. [qfactor()], [qinteger()],
+and [qdirichlet()] for the underlying transformations.
+}
+\author{
+Rob Carnell
+}
+\keyword{design}
diff --git a/tests/testthat/test-lhs_design.R b/tests/testthat/test-lhs_design.R
new file mode 100644
index 0000000..a868fae
--- /dev/null
+++ b/tests/testthat/test-lhs_design.R
@@ -0,0 +1,109 @@
+# Copyright 2026 Robert Carnell
+
+context("test-lhs_design")
+
+test_that("transform_lhs transforms each kind of margin", {
+ set.seed(1234)
+ # 5 single-column variables + a 3-column dirichlet = 8 LHS columns
+ X <- randomLHS(50, 8)
+ vars <- list(
+ x = list(dist = "qnorm", mean = 2, sd = 0.5),
+ grp = list(type = "factor", levels = c("a", "b", "c")),
+ ogrp = list(type = "factor", levels = c("lo", "mid", "hi"), ordered = TRUE),
+ count = list(type = "integer", min = 5L, max = 17L),
+ z = function(p) qexp(p, rate = 2),
+ comp = list(type = "dirichlet", alpha = c(2, 3, 4))
+ )
+ D <- transform_lhs(X, vars)
+
+ expect_s3_class(D, "data.frame")
+ expect_equal(nrow(D), 50)
+ # 5 single columns + dirichlet (3) - but comp consumes 3 columns and yields 3
+ expect_equal(ncol(D), 8)
+ expect_equal(names(D),
+ c("x", "grp", "ogrp", "count", "z", "comp1", "comp2", "comp3"))
+
+ # normal margin
+ expect_true(is.numeric(D$x))
+ # unordered factor
+ expect_true(is.factor(D$grp) && !is.ordered(D$grp))
+ expect_setequal(levels(D$grp), c("a", "b", "c"))
+ # ordered factor with levels in order
+ expect_true(is.ordered(D$ogrp))
+ expect_equal(levels(D$ogrp), c("lo", "mid", "hi"))
+ # integer in range
+ expect_true(all(D$count >= 5L & D$count <= 17L))
+ # function transform matches a direct call
+ expect_equal(D$z, qexp(X[, 5], rate = 2))
+ # dirichlet rows sum to one
+ expect_true(all(abs(rowSums(D[, c("comp1", "comp2", "comp3")]) - 1) < 1e-9))
+})
+
+test_that("lhs_design generates and transforms for each type", {
+ vars <- list(
+ a = list(dist = "qunif", min = 0, max = 10),
+ b = list(type = "integer", min = 1L, max = 6L)
+ )
+ for (ty in c("random", "improved", "maximin", "optimum", "genetic"))
+ {
+ set.seed(99)
+ D <- lhs_design(12, vars, type = ty)
+ expect_s3_class(D, "data.frame")
+ expect_equal(nrow(D), 12)
+ expect_equal(names(D), c("a", "b"))
+ expect_true(all(D$a >= 0 & D$a <= 10))
+ expect_true(all(D$b >= 1L & D$b <= 6L))
+ }
+
+ # oalhs may return a different (valid OA) number of rows
+ set.seed(99)
+ D <- lhs_design(9, vars, type = "oalhs")
+ expect_s3_class(D, "data.frame")
+ expect_equal(ncol(D), 2)
+})
+
+test_that("lhs_design passes extra arguments to the generator", {
+ set.seed(7)
+ vars <- list(a = function(p) p)
+ expect_silent(lhs_design(10, vars, type = "improved", dup = 3))
+ expect_silent(lhs_design(10, vars, type = "optimum", maxSweeps = 2, eps = 0.1))
+})
+
+test_that("lhs_design and transform_lhs validate their inputs", {
+ vars <- list(a = function(p) p, b = function(p) p)
+
+ # n validation
+ expect_error(lhs_design(0, vars))
+ expect_error(lhs_design(c(1, 2), vars))
+ expect_error(lhs_design(NA, vars))
+ expect_error(lhs_design(2.5, vars))
+
+ # variables validation
+ expect_error(lhs_design(10, list()))
+ expect_error(lhs_design(10, list(function(p) p))) # unnamed
+ expect_error(lhs_design(10, list(a = function(p) p, a = function(p) p)))
+ expect_error(lhs_design(10, list(a = "not a spec")))
+ expect_error(lhs_design(10, list(a = list(type = "factor")))) # no levels
+ expect_error(lhs_design(10, list(a = list(type = "integer")))) # no min/max
+ expect_error(lhs_design(10, list(a = list(type = "dirichlet")))) # no alpha
+
+ # transform_lhs column count must match
+ X <- randomLHS(10, 3)
+ expect_error(transform_lhs(X, vars)) # needs 2 columns, given 3
+ expect_error(transform_lhs(as.data.frame(X), vars)) # not a matrix
+ Y <- randomLHS(10, 2)
+ Y[1, 1] <- NA
+ expect_error(transform_lhs(Y, vars))
+ Z <- randomLHS(10, 2)
+ Z[1, 1] <- 2
+ expect_error(transform_lhs(Z, vars))
+})
+
+test_that("transform_lhs works on augmented designs", {
+ set.seed(321)
+ X <- augmentLHS(randomLHS(6, 2), 4)
+ vars <- list(a = list(dist = "qnorm"), b = list(type = "integer", min = 0L, max = 9L))
+ D <- transform_lhs(X, vars)
+ expect_equal(nrow(D), 10)
+ expect_equal(names(D), c("a", "b"))
+})