logLikGrad: Log-Likelihood Function

View source: R/logLikGrad.R

logLikGradR Documentation

Log-Likelihood Function

Description

Calculates the gradient of the negative “concentrated” log-likelihood.

Usage

logLikGrad(param, object, ...)

## Default S3 method:
logLikGrad(param, object, x, y, dx = NULL, dy = NULL, 
	covtype = c("matern5_2", "matern3_2", "gaussian"),
	tolerance = NULL, envir, ...)
## S3 method for class 'gekm'
logLikGrad(param, object, ...)

Arguments

param

a numeric vector corresponding to the values of the correlation parameters at which the negative “concentrated” log-likelihood is to be evaluated.

object

a numeric matrix containing the or an object of class gekm.

x

a model.matrix.

y

a vector of response values.

dx

the derivatives of the model matrix x, see derivModelMatrix for details. Default is NULL.

dy

the vector of derivatives of the response y. Default is NULL.

covtype

a character to specify the correlation structure to be used. One of matern5_2, matern3_2 or gaussian. Default is matern5_2, see blockCor for details.

tolerance

a tolerance for the conditional number of the correlation matrix, see blockChol for details. Default is NULL, i.e. no regularization is applied.

envir

an environment with auxiliary variables obtained from a corresponding call to logLikFun.

...

arguments to be passed to the default method.

Value

The value of the negative “concentrated” log-likelihood at param.

Author(s)

Carmen van Meegen

References

Park, J.-S. and Beak, J. (2001). Efficient Computation of Maximum Likelihood Estimators in a Spatial Linear Model with Power Exponential Covariogram. Computers & Geosciences, 27(1):1–7. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1016/S0098-3004(00)00016-9")}.

Rasmussen, C. E. and Williams, C. K. I. (2006). Gaussian Processes for Machine Learning. The MIT Press. https://gaussianprocess.org/gpml/.

Santner, T. J., Williams, B. J., and Notz, W. I. (2018). The Design and Analysis of Computer Experiments. 2nd edition. Springer-Verlag.

Zimmermann, R. (2015). On the Condition Number Anomaly of Gaussian Correlation Matrices. Linear Algebra and its Applications, 466:512–526. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1016/j.laa.2014.10.038")}.

See Also

logLikFun for computing the value of the negative “concentrated” log-likelihood.

gekm for fitting a (gradient-enhanced) Kriging model.

Examples

## 2-dimensional example

# Generate coordinates and calculate slopes
x1 <- seq(-1.75, 1.75, length = 3)
x2 <- seq(-0.75, 0.75, length = 3)
X <- expand.grid(x1 = x1, x2 = x2)
y <- camel6(X)
dy <- camel6Grad(X)
dat <- data.frame(X, y)
deri <- data.frame(dy)

# Fit (gradient-enhanced) Kriging model
km.2d <- gekm(y ~ 1, data = dat, covtype = "gaussian", optimizer = "L-BFGS-B")
gekm.2d <- gekm(y ~ 1, data = dat, deriv = deri, covtype = "gaussian", optimizer = "L-BFGS-B")

# Compute negative 'concentrated' log-likelihood values
n.grid <- 20
theta1.grid <- seq(0.5, 4, length = n.grid)
theta2.grid <- seq(0.5, 2, length = n.grid)
params <- expand.grid(theta1 = theta1.grid, theta2 = theta2.grid)

logLik.km.2d <- apply(params, 1, logLikFun, km.2d)
logLik.gekm.2d <- apply(params, 1, logLikFun, gekm.2d)

logLikGrad.km.2d <- t(apply(params, 1, logLikGrad, km.2d))
logLikGrad.gekm.2d <- t(apply(params, 1, logLikGrad, gekm.2d))

# Plot negative 'concentrated' log-likelihood
par(mfrow = c(1, 2), oma = c(3.6, 3.5, 1.5, 0.2), mar = c(0, 0, 1.5, 0))
contour(theta1.grid, theta2.grid, matrix(logLik.km.2d, nrow = n.grid, ncol = n.grid), 
		nlevels = 50, main = "Kriging")
vectorfield(params, logLikGrad.km.2d, col = 4, lwd = 2, length = 0.1)
points(km.2d$theta[1], km.2d$theta[2], col = "red", pch = 16)

contour(theta1.grid, theta2.grid, matrix(logLik.gekm.2d, nrow = n.grid, ncol = n.grid), 
		nlevels = 50, main = "GEK", yaxt = "n")
points(gekm.2d$theta[1], gekm.2d$theta[2], col = "red", pch = 16)
vectorfield(params, logLikGrad.gekm.2d, col = 4, lwd = 2, length = 0.1)

title(main = "Negative 'concentrated' log-likelihood", outer = TRUE)
mtext(side = 1, outer = TRUE, line = 2.5, expression(theta[1]))
mtext(side = 2, outer = TRUE, line = 2.5, expression(theta[2]))

gek documentation built on Jan. 31, 2026, 1:07 a.m.

Related to logLikGrad in gek...