0.6.2-dev0
FORCES
FORtran lib for Comp. Env. Sys.
Loading...
Searching...
No Matches
mo_anneal::anneal Interface Reference

Optimize cost function with simulated annealing. More...

Public Member Functions

real(dp) function, dimension(size(para, 1)) anneal_dp (eval, cost, para, prange, prange_func, temp, dt, nitermax, len, nst, eps, acc, seeds, printflag, maskpara, weight, changeparamode, reflectionflag, pertubflexflag, maxit, undef_funcval, tmp_file, funcbest, history)
 

Detailed Description

Optimize cost function with simulated annealing.

Optimizes a user provided cost function using the Simulated Annealing strategy.

Example

User defined function cost_dp which calculates the cost function value for a parameter set (the interface given below has to be used for this function!).

para = (/ 1.0_dp , 2.0_dp /)
prange(1,:) = (/ 0.0_dp, 10.0_dp /)
prange(2,:) = (/ 0.0_dp, 50.0_dp /)
parabest = anneal(cost_dp, para, prange)

See also test folder for a detailed example, "test/test_mo_anneal".

Literature

  1. S. Kirkpatrick, C. D. Gelatt, and M. P. Vecchi. Optimization by simulated annealing. Science, 220:671-680, 1983.
  2. N. Metropolis, A. W. Rosenbluth, M. N. Rosenbluth, A. H. Teller, and E. Teller. Equation of state calculations by fast computing machines. J. Chem. Phys., 21:1087-1092, June 1953.
  3. B. A. Tolson, and C. A. Shoemaker. Dynamically dimensioned search algorithm for computationally efficient watershed model calibration. WRR, 43(1), W01413, 2007.
Parameters
[in]INTERFACE :: evalInterface calculating the eval function at a given point.
[in]INTERFACE :: costInterface calculating the cost function at a given point.
[in]real(dp), dimension(:) :: paraInitial parameter set.
[in]real(dp), dimension(size(para),2), optional :: prangeLower and upper bound per parameter.
[in]interface, optional :: prange_funcInterface calculating the feasible range for a parameter at a certain point, if ranges are variable.
[in]real(dp), optional :: tempInitial temperature.
DEFAULT: Get_Temperature.
[in]real(dp), optional :: DTGeometrical decreement of temperature.
0.7<DT<0.999.
DEFAULT: 0.9_dp.
[in]integer(i4), optional :: nITERmaxMaximal number of iterations.
Will be increased by 10% if stopping criteria of acceptance ratio or epsilon decreement of cost function is not fullfilled.
DEFAULT: 1000_i4.
[in]integer(i4), optional :: LENLength of Markov Chain.
DEFAULT: MAX(250_i4, size(para,1)).
[in]integer(i4), optional :: nSTNumber of consecutive LEN steps.
DEFAULT: 5_i4.
[in]real(dp), optional :: epsStopping criteria of epsilon decreement of cost function
DEFAULT: 0.0001_dp
[in]real(dp), optional :: accStopping criteria for Acceptance Ratio
acc <= 0.1_dp
DEFAULT: 0.1_dp
[in]integer(i4/i8), dimension(3), optional :: seedsSeeds of random numbers used for random parameter set generation
DEFAULT: dependent on current time
[in]logical, optional :: printflagIf .true. detailed command line output is written
DEFAULT: .false.
[in]logical, dimension(size(para)), optional :: maskparamaskpara(i) = .true. --> parameter is optimized
maskpara(i) = .false. --> parameter is discarded from optimiztaion
DEFAULT: .true.
[in]real(dp), dimension(size(para)), optional :: weightvector of weights per parameter:
gives the frequency of parameter to be chosen for optimization (will be scaled to a CDF internally)
eg. [1,2,1] --> parameter 2 is chosen twice as often as parameter 1 and 2
DEFAULT: weight = 1.0_dp
[in]integer(i4), optional :: changeParaModewhich and how many param.are changed in one step
1 = one parameter
2 = all parameter
3 = neighborhood parameter
DEFAULT: 1_i4
[in]logical, optional :: reflectionFlagif new parameter values are Gaussian distributed and reflected (.true.) or uniform in range (.false.)
DEFAULT: .false.
[in]logical, optional :: pertubFlexFlagif pertubation of Gaussian distributed parameter values is constant at 0.2 (.false.) or depends on dR (.true.)
DEFAULT: .true.
[in]logical, optional :: maxitmaximizing (.true.) or minimizing (.false.) a function
DEFAULT: .false. (minimization)
[in]real(dp), optional :: undef_funcvalobjective function value defining invalid model output, e.g. -999.0_dp
[in]character(len=*) , optional :: tmp_filefile with temporal output
[out]real(dp), optional :: funcbestminimized value of cost function
[out]real(dp), dimension(:,:), allocatable, optional :: historyreturns a vector of achieved objective after ith model evaluation
Return values
real(dp) :: parabest(size(para))Parameter set minimizing the cost function.
Note
  • Either fixed parameter range (prange) OR flexible parameter range (function interface prange_func) has to be given in calling sequence.
  • Only double precision version available.
  • If single precision is needed not only dp has to be replaced by sp but also i8 of save_state (random number variables) has to be replaced by i4.
  • ParaChangeMode > 1 is not applied in GetTemperature.
  • For Temperature estimation always only one single parameter is changed (ParaChangeMode=1) which should give theoretically always the best estimate.
  • cost_func and prange_func are user defined functions. See interface definition.
Author
Luis Samaniego
Date
Jan 2000
Mar 2003
  • Re-heating
Author
Juliane Mai
Date
Mar 2012
  • Modular version
May 2012
  • sp version
  • documentation

Definition at line 186 of file mo_anneal.f90.

Member Function/Subroutine Documentation

◆ anneal_dp()

real(dp) function, dimension(size(para, 1)) mo_anneal::anneal::anneal_dp ( procedure(eval_interface), intent(in), pointer  eval,
procedure(objective_interface), intent(in), pointer  cost,
real(dp), dimension(:), intent(in)  para,
real(dp), dimension(size(para, 1), 2), intent(in), optional  prange,
external subroutine(real(dp), dimension(:), intent(in) paraset, integer(i4), intent(in) ipar, real(dp), dimension(2), intent(out) rangepar), optional  prange_func,
real(dp), intent(in), optional  temp,
real(dp), intent(in), optional  dt,
integer(i4), intent(in), optional  nitermax,
integer(i4), intent(in), optional  len,
integer(i4), intent(in), optional  nst,
real(dp), intent(in), optional  eps,
real(dp), intent(in), optional  acc,
integer(i8), dimension(3), intent(in), optional  seeds,
logical, intent(in), optional  printflag,
logical, dimension(size(para, 1)), intent(in), optional  maskpara,
real(dp), dimension(size(para, 1)), intent(in), optional  weight,
integer(i4), intent(in), optional  changeparamode,
logical, intent(in), optional  reflectionflag,
logical, intent(in), optional  pertubflexflag,
logical, intent(in), optional  maxit,
real(dp), intent(in), optional  undef_funcval,
character(len = *), intent(in), optional  tmp_file,
real(dp), intent(out), optional  funcbest,
real(dp), dimension(:, :), intent(out), optional, allocatable  history 
)
Parameters
[in]parainitial parameter
[in]prangelower and upper limit per parameter
[in]tempstarting temperature (DEFAULT: Get_Temperature)
[in]dtgeometrical decreement, 0.7<DT<0.999 (DEFAULT: 0.9)
[in]nitermaxmaximal number of iterations (DEFAULT: 1000)
[in]lenLength of Markov Chain, DEFAULT: max(250, size(para,1))
[in]nstNumber of consecutive LEN steps! (DEFAULT: 5)
[in]epsepsilon decreement of cost function (DEFAULT: 0.01)
[in]accAcceptance Ratio, <0.1 stopping criteria (DEFAULT: 0.1)
[in]seedsSeeds of random numbers (DEFAULT: Get_timeseed)
[in]printflagIf command line output is written (.true.) (DEFAULT: .false.)
[in]maskparatrue if parameter will be optimized false if parameter is discarded in optimization (DEFAULT: .true.)
[in]weightvector of weights per parameter gives the frequency of parameter to be chosen for optimization (DEFAULT: uniform)
[in]changeparamodewhich and how many param. are changed in a step 1 = one parameter 2 = all parameter 3 = neighborhood parameter (DEFAULT: 1_i4)
[in]reflectionflagif new parameter values are selected normal distributed and reflected (.true.) or uniform in range (.false.) (DEFAULT: .false.)
[in]pertubflexflagif pertubation of normal distributed parameter values is constant 0.2 (.false.) or depends on dR (.true.) (DEFAULT: .true.)
[in]maxitMaximization or minimization of function maximization = .true., minimization = .false. (DEFAULT: .false.)
[in]undef_funcvalobjective function value occuring if parameter set leads to invalid model results, e.g. -9999.0_dp! (DEFAULT: not present)
[in]tmp_filefile for temporal output
[out]funcbestminimized value of cost function (DEFAULT: not present)
[out]historyreturns a vector of achieved objective! after ith model evaluation (DEFAULT: not present)
Returns
parameter set minimizing objective

Definition at line 289 of file mo_anneal.f90.


The documentation for this interface was generated from the following file: