ForK Library
Functions/Subroutines

krigingextremefuncpredict.f90 File Reference

Subroutine used to determine the minimum or maximum value from a Kriging response surface (requires buildkriging to be called first)
Only Works for Ordinary Kriging (Constant Mean function) More...

Go to the source code of this file.

Functions/Subroutines

subroutine krigingextremefuncpredict (ndim, ntot, X, stot, H, beta, V, hyper, Xm, Ym, covarflagi, optflag, flag, Lb, Ub)
 This subroutine predicts the global minimum or maximum function values from the Kriging surface (constructed using buildkriging). Using an existing Kriging surface, the minimum or maximum value from the model is determined using either simplexsearch or patternsearch (specified with the optflagi). Within these optimization subroutines, the code krigingfuncpredict is invoked to predict the mean function value of the Kriging surface at a particular location. The optimization algorithms by design determine the minimum value of a function. To determine the maximum value (specified by flag=1), the objective is multiplied by \(-1\) and minimization is performed.

This optimization can only be performed on an Ordinary Kriging surface, meaning that the mean function is a constant value that is fitted during the construction of the mean function. From the universal Kriging framework used in buildkriging, ordinary Kriging can be recovered by using a \(p=0\) regression. This zeroth order regression requires a single term (meaning \(stot=1\)) and all the elements of the collocation matrix (H) are equal to one.

Detailed Description

Subroutine used to determine the minimum or maximum value from a Kriging response surface (requires buildkriging to be called first)
Only Works for Ordinary Kriging (Constant Mean function)

Definition in file krigingextremefuncpredict.f90.


Function Documentation

subroutine krigingextremefuncpredict ( integer,intent(in)  ndim,
integer,intent(in)  ntot,
real(8),dimension(ndim,ntot),intent(in)  X,
integer,intent(in)  stot,
real(8),dimension(stot,ntot),intent(in)  H,
real(8),dimension(stot),intent(in)  Beta,
real(8),dimension(ntot),intent(in)  V,
real(8),dimension(ndim+2),intent(in)  hyper,
real(8),dimension(ndim),intent(out)  Xm,
real(8),intent(out)  Ym,
integer,intent(in)  covarflagi,
integer,intent(in)  optflag,
integer,intent(in)  flag,
real(8),dimension(ndim),intent(in)  Lb,
real(8),dimension(ndim),intent(in)  Ub 
)

This subroutine predicts the global minimum or maximum function values from the Kriging surface (constructed using buildkriging). Using an existing Kriging surface, the minimum or maximum value from the model is determined using either simplexsearch or patternsearch (specified with the optflagi). Within these optimization subroutines, the code krigingfuncpredict is invoked to predict the mean function value of the Kriging surface at a particular location. The optimization algorithms by design determine the minimum value of a function. To determine the maximum value (specified by flag=1), the objective is multiplied by \(-1\) and minimization is performed.

This optimization can only be performed on an Ordinary Kriging surface, meaning that the mean function is a constant value that is fitted during the construction of the mean function. From the universal Kriging framework used in buildkriging, ordinary Kriging can be recovered by using a \(p=0\) regression. This zeroth order regression requires a single term (meaning \(stot=1\)) and all the elements of the collocation matrix (H) are equal to one.

Author:
Brian Lockwood Department of Mechanical Engineering University of Wyoming
Date:
May 2, 2012
Parameters:
in)ndim : The dimension of the problem
in)ntot : The number of Training points
in)X : The location of the training points (size=[ndimxntot])
in)stot : Number of Terms in the regression
in)H: The collocation matrix for the regression (size=[stotxntot])
in) beta: Regression coefficients based on the optimal estimate for the Kriging model (size=[stot])
Supplied by buildkriging subroutine
in)V: Processed Training Data (size=[ntot])
Supplied by buildkriging subroutine
in)hyper: Hyperparameters for the Kriging Model (size=[ndim+2])
Supplied by buildkriging subroutine
inout) Xm : Location of the extreme point (size=[ndim])
Starting point for the optimization is the input value for Xm. Final value overwrites this initial value.
out) Ym : Minimum or Maximum value
in)covarflagi: Flag to govern which covariance function is used

  • covarflag==0 Uses Matern function with \(\nu=1/2\)
  • covarflag==1 Uses Matern function with \(\nu=3/2\)
  • covarflag==2 Uses Matern function with \(\nu=5/2\)


The parameter \(\nu\) governs the smoothness and differentiability of the covariance function. For function only Kriging, all three options are available.
Good rule of thumb is to use the least smooth ( \(\nu=1/2\)) unless there is a reason to assume more smoothness in the data
For small numbers of training points, higher \(\nu\) can improve accuracy

in)optflagi: Flag to govern the optimization algorithm used to determine minimum or maximum value

  • optflag=0 Simplex Search (Nelder-Mead Method) to determine extreme value
    Local optimization technique using simplex elements to determine search direction and line search to determine step sizes
    Fastest method but inherently sequential and may stall on local optimum
  • optflag=1 Pattern Search used to determine extreme value
    Global optimization method with fixed search directions (forward and backward in each direction)
    Each iteration requires \(2*ndim\) function evaluations; however, these may be performed in parallel using openmp (Each thread performs own likelihood evaluation so OMP_STACK_SIZE must be set large enough)
in) flag Determines whether minimum or maximum value is found

  • flag = 0 Subroutine finds the minimum value
  • flag = 1 Subroutine finds the maximum value
in) Lb Lower Bound for each variable (size = [ndim])
Should be the minimum possible value of Xm in each dimension
in) Ub Upper Bound for each variable (size = [ndim])
Should be the maximum possible value of Xm in each dimension

Definition at line 57 of file krigingextremefuncpredict.f90.

References patternsearch(), and simplexsearch().

 All Classes Namespaces Files Functions Variables