0.6.2-dev0
FORCES
FORtran lib for Comp. Env. Sys.
Loading...
Searching...
No Matches
mo_nelmin Module Reference

Nelder-Mead algorithm. More...

Data Types

interface  nelminrange
 Minimizes a user-specified function using the Nelder-Mead algorithm. More...
 

Functions/Subroutines

real(dp) function, dimension(size(pstart)), public nelmin (func, pstart, varmin, step, konvge, maxeval, funcmin, neval, numrestart, ierror, history)
 Minimizes a user-specified function using the Nelder-Mead algorithm.
 
real(dp) function, dimension(size(pstart)), public nelminxy (func, pstart, xx, yy, varmin, step, konvge, maxeval, funcmin, neval, numrestart, ierror, history)
 Minimizes a user-specified function using the Nelder-Mead algorithm.
 
real(dp) function, dimension(size(pstart)) nelminrange_dp (func, pstart, prange, varmin, step, konvge, maxeval, funcmin, neval, numrestart, ierror, history)
 
real(sp) function, dimension(size(pstart)) nelminrange_sp (func, pstart, prange, varmin, step, konvge, maxeval, funcmin, neval, numrestart, ierror, history)
 

Detailed Description

Nelder-Mead algorithm.

This module provides NELMIN, which minimizes a function using the Nelder-Mead algorithm with the Applied Statistics algorithms No. 047. Original FORTRAN77 version by R ONeill. FORTRAN90 version by John Burkardt.

Author
Matthias Cuntz
Date
2012

Function/Subroutine Documentation

◆ nelmin()

real(dp) function, dimension(size(pstart)), public mo_nelmin::nelmin ( external real(dp) function(real(dp), dimension(:), intent(in) pp)  func,
real(dp), dimension(:), intent(in)  pstart,
real(dp), intent(in), optional  varmin,
real(dp), dimension(:), intent(in), optional  step,
integer(i4), intent(in), optional  konvge,
integer(i4), intent(in), optional  maxeval,
real(dp), intent(out), optional  funcmin,
integer(i4), intent(out), optional  neval,
integer(i4), intent(out), optional  numrestart,
integer(i4), intent(out), optional  ierror,
real(dp), dimension(:), intent(out), optional, allocatable  history 
)

Minimizes a user-specified function using the Nelder-Mead algorithm.

Simplex function minimisation procedure due to Nelder and Mead (1965), as implemented by O'Neill(1971, Appl.Statist. 20, 338-45), with subsequent comments by Chambers+Ertel(1974, 23, 250-1), Benyon(1976, 25, 97) and Hill(1978, 27, 380-2)

The function to be minimized is the first argument of nelmin and must be defined as

FUNCTION func(p)
USE mo_kind, ONLY: dp
IMPLICIT NONE
REAL(dp), DIMENSION(:), INTENT(IN) :: p
REAL(dp) :: func
END FUNCTION func
Define number representations.
Definition mo_kind.F90:17
integer, parameter dp
Double Precision Real Kind.
Definition mo_kind.F90:46

and for nelminxy

FUNCTION func(p,x,y)
USE mo_kind, ONLY: dp
IMPLICIT NONE
REAL(dp), DIMENSION(:), INTENT(IN) :: p
REAL(dp), DIMENSION(:), INTENT(IN) :: x
REAL(dp), DIMENSION(:), INTENT(IN) :: y
REAL(dp) :: func
END FUNCTION func

This routine does not include a termination test using the fitting of a quadratic surface.

Calling sequence

pmin = nelmin(func, pstart, funcmin, varmin, step, konvge, maxeval, &
neval, numrestart, ierror)
pmin = nelminxy(func, pstart, xx, yy, funcmin, varmin, step, konvge, maxeval, &
neval, numrestart, ierror)
pmin = nelmin(func, pstart, prange, funcmin, varmin, step, konvge, maxeval, &
neval, numrestart, ierror)

Original FORTRAN77 version by R ONeill. FORTRAN90 version by John Burkardt.

Example

pstart(1:n) = (/ -1.2_dp, 1.0_dp /)
step(1:n) = (/ 1.0_dp, 1.0_dp /)
pmin = nelmin(rosenbrock, pstart, step=step, ierror=ierror)

-> see also example in test directory

Literature

  1. Simplex function minimisation procedure due to Nelder and Mead (1965), as implemented by O'Neill(1971, Appl.Statist. 20, 338-45), with subsequent comments by Chambers+Ertel(1974, 23, 250-1), Benyon(1976, 25, 97) and Hill(1978, 27, 380-2)
  2. John Nelder, Roger Mead, A simplex method for function minimization, Computer Journal, Volume 7, 1965, pages 308-313.
  3. R ONeill, Algorithm AS 47: Function Minimization Using a Simplex Procedure, Applied Statistics, Volume 20, Number 3, 1971, pages 338-345.
Parameters
[in]real(dp) :: func(p,xx,yy)Function on which to search the minimum
[in]real(dp) :: pstart(:)Starting point for the iteration.
[in]real(dp) :: xx(ONLY nelminxy) First values to pass as function arguments
[in]real(dp) :: yy(ONLY nelminxy) Second values to pass as function arguments
[in]real(dp) :: prange(:,2)(ONLY nelminrange) Range of parameters (upper and lower bound).
[in]real(dp), optional :: varminthe terminating limit for the variance of the function values. varmin>0 is required.
[in]real(dp), optional :: step(:)determines the size and shape of the initial simplex. The relative magnitudes of its elements should reflect the units of the variables. size(step)=size(start)
[in]integer(i4), optional :: konvgethe convergence check is carried out every konvge iterations. konvge>0 is required.
[in]integer(i4), optional :: maxevalthe maximum number of function evaluations. default: 1000
[out]real(dp), optional :: funcminthe minimum value of the function.
[out]integer(i4), optional :: nevalthe number of function evaluations used.
[out]integer(i4), optional :: numrestartthe number of restarts.
[out]integer(i4), optional :: ierrorerror indicator. 0: no errors detected. 1: varmin or konvge have an illegal value. 2: terminated because maxeval was exceeded without convergence.
[out]real(dp), optional, allocatable :: history(:)the history of best function values, history(neval)=funcmin
Returns
real(sp/dp) :: nelmin(size(start)) — coordinates of the point which is estimated to minimize the function.
Author
Matthias Cuntz
Date
Jul 2012
  • i4, dp, intent
  • function, optional
  • nelminxy
Author
Juliane Mai
Date
Aug 2012
  • nelminrange
Dec 2012
  • history output

Definition at line 136 of file mo_nelmin.f90.

References mo_kind::dp, and nelmin().

Referenced by nelmin().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ nelminrange_dp()

real(dp) function, dimension(size(pstart)) mo_nelmin::nelminrange_dp ( external real(dp) function(real(dp), dimension(:), intent(in) pp)  func,
real(dp), dimension(:), intent(in)  pstart,
real(dp), dimension(:,:), intent(in)  prange,
real(dp), intent(in), optional  varmin,
real(dp), dimension(:), intent(in), optional  step,
integer(i4), intent(in), optional  konvge,
integer(i4), intent(in), optional  maxeval,
real(dp), intent(out), optional  funcmin,
integer(i4), intent(out), optional  neval,
integer(i4), intent(out), optional  numrestart,
integer(i4), intent(out), optional  ierror,
real(dp), dimension(:), intent(out), optional, allocatable  history 
)
private

Definition at line 817 of file mo_nelmin.f90.

◆ nelminrange_sp()

real(sp) function, dimension(size(pstart)) mo_nelmin::nelminrange_sp ( external real(sp) function(real(sp), dimension(:), intent(in) pp)  func,
real(sp), dimension(:), intent(in)  pstart,
real(sp), dimension(:,:), intent(in)  prange,
real(sp), intent(in), optional  varmin,
real(sp), dimension(:), intent(in), optional  step,
integer(i4), intent(in), optional  konvge,
integer(i4), intent(in), optional  maxeval,
real(sp), intent(out), optional  funcmin,
integer(i4), intent(out), optional  neval,
integer(i4), intent(out), optional  numrestart,
integer(i4), intent(out), optional  ierror,
real(sp), dimension(:), intent(out), optional, allocatable  history 
)
private

Definition at line 1174 of file mo_nelmin.f90.

◆ nelminxy()

real(dp) function, dimension(size(pstart)), public mo_nelmin::nelminxy ( external real(dp) function(real(dp), dimension(:), intent(in) pp, real(dp), dimension(:), intent(in) xxx, real(dp), dimension(:), intent(in) yyy)  func,
real(dp), dimension(:), intent(in)  pstart,
real(dp), dimension(:), intent(in)  xx,
real(dp), dimension(:), intent(in)  yy,
real(dp), intent(in), optional  varmin,
real(dp), dimension(:), intent(in), optional  step,
integer(i4), intent(in), optional  konvge,
integer(i4), intent(in), optional  maxeval,
real(dp), intent(out), optional  funcmin,
integer(i4), intent(out), optional  neval,
integer(i4), intent(out), optional  numrestart,
integer(i4), intent(out), optional  ierror,
real(dp), dimension(:), intent(out), optional, allocatable  history 
)

Minimizes a user-specified function using the Nelder-Mead algorithm.

Simplex function minimisation procedure due to Nelder and Mead (1965), as implemented by O'Neill(1971, Appl.Statist. 20, 338-45), with subsequent comments by Chambers+Ertel(1974, 23, 250-1), Benyon(1976, 25, 97) and Hill(1978, 27, 380-2)

The function to be minimized is the first argument of nelmin and must be defined as

FUNCTION func(p)
USE mo_kind, ONLY: dp
IMPLICIT NONE
REAL(dp), DIMENSION(:), INTENT(IN) :: p
REAL(dp) :: func
END FUNCTION func

and for nelminxy

FUNCTION func(p,x,y)
USE mo_kind, ONLY: dp
IMPLICIT NONE
REAL(dp), DIMENSION(:), INTENT(IN) :: p
REAL(dp), DIMENSION(:), INTENT(IN) :: x
REAL(dp), DIMENSION(:), INTENT(IN) :: y
REAL(dp) :: func
END FUNCTION func

This routine does not include a termination test using the fitting of a quadratic surface.

Calling sequence

pmin = nelmin(func, pstart, funcmin, varmin, step, konvge, maxeval, &
neval, numrestart, ierror)
pmin = nelminxy(func, pstart, xx, yy, funcmin, varmin, step, konvge, maxeval, &
neval, numrestart, ierror)
pmin = nelmin(func, pstart, prange, funcmin, varmin, step, konvge, maxeval, &
neval, numrestart, ierror)

Original FORTRAN77 version by R ONeill. FORTRAN90 version by John Burkardt.

Example

pstart(1:n) = (/ -1.2_dp, 1.0_dp /)
step(1:n) = (/ 1.0_dp, 1.0_dp /)
pmin = nelmin(rosenbrock, pstart, step=step, ierror=ierror)

-> see also example in test directory

Literature

  1. Simplex function minimisation procedure due to Nelder and Mead (1965), as implemented by O'Neill(1971, Appl.Statist. 20, 338-45), with subsequent comments by Chambers+Ertel(1974, 23, 250-1), Benyon(1976, 25, 97) and Hill(1978, 27, 380-2)
  2. John Nelder, Roger Mead, A simplex method for function minimization, Computer Journal, Volume 7, 1965, pages 308-313.
  3. R ONeill, Algorithm AS 47: Function Minimization Using a Simplex Procedure, Applied Statistics, Volume 20, Number 3, 1971, pages 338-345.
Parameters
[in]real(dp) :: func(p,xx,yy)Function on which to search the minimum
[in]real(dp) :: pstart(:)Starting point for the iteration.
[in]real(dp) :: xx(ONLY nelminxy) First values to pass as function arguments
[in]real(dp) :: yy(ONLY nelminxy) Second values to pass as function arguments
[in]real(dp) :: prange(:,2)(ONLY nelminrange) Range of parameters (upper and lower bound).
[in]real(dp), optional :: varminthe terminating limit for the variance of the function values. varmin>0 is required.
[in]real(dp), optional :: step(:)determines the size and shape of the initial simplex. The relative magnitudes of its elements should reflect the units of the variables. size(step)=size(start)
[in]integer(i4), optional :: konvgethe convergence check is carried out every konvge iterations. konvge>0 is required.
[in]integer(i4), optional :: maxevalthe maximum number of function evaluations. default: 1000
[out]real(dp), optional :: funcminthe minimum value of the function.
[out]integer(i4), optional :: nevalthe number of function evaluations used.
[out]integer(i4), optional :: numrestartthe number of restarts.
[out]integer(i4), optional :: ierrorerror indicator. 0: no errors detected. 1: varmin or konvge have an illegal value. 2: terminated because maxeval was exceeded without convergence.
[out]real(dp), optional, allocatable :: history(:)the history of best function values, history(neval)=funcmin
Returns
real(sp/dp) :: nelmin(size(start)) — coordinates of the point which is estimated to minimize the function.
Author
Matthias Cuntz
Date
Jul 2012
  • i4, dp, intent
  • function, optional
  • nelminxy
Author
Juliane Mai
Date
Aug 2012
  • nelminrange
Dec 2012
  • history output

Definition at line 475 of file mo_nelmin.f90.

References mo_kind::dp, nelminxy(), and mo_kind::sp.

Referenced by nelminxy().

Here is the call graph for this function:
Here is the caller graph for this function: