FORCES
FORtran lib for Comp. Env. Sys.
|
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) |
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.
COPYING
and COPYING.LESSER
provided with this software. The complete GNU license text can also be found at http://www.gnu.org/licenses/. 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
and for nelminxy
This routine does not include a termination test using the fitting of a quadratic surface.
Calling sequence
Original FORTRAN77 version by R ONeill. FORTRAN90 version by John Burkardt.
Example
-> see also example in test directory
Literature
[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 :: varmin | the 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 :: konvge | the convergence check is carried out every konvge iterations. konvge>0 is required. |
[in] | integer(i4), optional :: maxeval | the maximum number of function evaluations. default: 1000 |
[out] | real(dp), optional :: funcmin | the minimum value of the function. |
[out] | integer(i4), optional :: neval | the number of function evaluations used. |
[out] | integer(i4), optional :: numrestart | the number of restarts. |
[out] | integer(i4), optional :: ierror | error 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 |
Definition at line 136 of file mo_nelmin.f90.
References mo_kind::dp, and nelmin().
Referenced by nelmin().
|
private |
Definition at line 817 of file mo_nelmin.f90.
|
private |
Definition at line 1174 of file mo_nelmin.f90.
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
and for nelminxy
This routine does not include a termination test using the fitting of a quadratic surface.
Calling sequence
Original FORTRAN77 version by R ONeill. FORTRAN90 version by John Burkardt.
Example
-> see also example in test directory
Literature
[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 :: varmin | the 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 :: konvge | the convergence check is carried out every konvge iterations. konvge>0 is required. |
[in] | integer(i4), optional :: maxeval | the maximum number of function evaluations. default: 1000 |
[out] | real(dp), optional :: funcmin | the minimum value of the function. |
[out] | integer(i4), optional :: neval | the number of function evaluations used. |
[out] | integer(i4), optional :: numrestart | the number of restarts. |
[out] | integer(i4), optional :: ierror | error 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 |
Definition at line 475 of file mo_nelmin.f90.
References mo_kind::dp, nelminxy(), and mo_kind::sp.
Referenced by nelminxy().