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

Shuffled Complex Evolution optimization algorithm. More...

Functions/Subroutines

real(dp) function, dimension(size(pini, 1)), public sce (eval, functn, pini, prange, mymaxn, mymaxit, mykstop, mypcento, mypeps, myseed, myngs, mynpg, mynps, mynspl, mymings, myiniflg, myprint, mymask, myalpha, mybeta, tmp_file, popul_file, popul_file_append, parallel, restart, restart_file, bestf, neval, history)
 Shuffled Complex Evolution (SCE) algorithm for global optimization.
 
subroutine parstt (x, bound, peps, mask, xnstd, gnrng, ipcnvg)
 
subroutine comp (ngs2, npg, a, af, b, bf)
 
subroutine sort_matrix (rb, ra)
 
subroutine chkcst (x, bl, bu, mask, ibound)
 
subroutine getpnt (idist, bl, bu, std, xi, mask, save_state, x)
 
subroutine cce (s, sf, bl, bu, maskpara, xnstd, icall, maxn, maxit, save_state_gauss, functn, eval, alpha, beta, history, idot ifdef mpi
 

Detailed Description

Shuffled Complex Evolution optimization algorithm.

Optimization algorithm using Shuffled Complex Evolution strategy. Original version 2.1 of Qingyun Duan (1992) rewritten in Fortran 90.

Authors
Juliane Mai
Date
Feb 2013

Function/Subroutine Documentation

◆ cce()

subroutine mo_sce::cce ( real(dp), dimension(:, :), intent(inout)  s,
real(dp), dimension(:), intent(inout)  sf,
real(dp), dimension(:), intent(in)  bl,
real(dp), dimension(:), intent(in)  bu,
logical, dimension(:), intent(in)  maskpara,
real(dp), dimension(:), intent(in)  xnstd,
integer(i8), intent(inout)  icall,
integer(i8), intent(in)  maxn,
logical, intent(in)  maxit,
integer(i8), dimension(n_save_state), intent(inout)  save_state_gauss,
procedure(objective_interface), intent(in), pointer  functn,
procedure(eval_interface), intent(in), pointer  eval,
real(dp), intent(in)  alpha,
real(dp), intent(in)  beta,
real(dp), dimension(maxn), intent(inout)  history,
logical, intent(in)  idot,
  ifdef,
  mpi 
)
private

Definition at line 1830 of file mo_sce.F90.

◆ chkcst()

subroutine mo_sce::chkcst ( real(dp), dimension(:), intent(in)  x,
real(dp), dimension(:), intent(in)  bl,
real(dp), dimension(:), intent(in)  bu,
logical, dimension(:), intent(in)  mask,
integer(i4), intent(out)  ibound 
)
private

Definition at line 1735 of file mo_sce.F90.

◆ comp()

subroutine mo_sce::comp ( integer(i4), intent(in)  ngs2,
integer(i4), intent(in)  npg,
real(dp), dimension(:, :), intent(inout)  a,
real(dp), dimension(:), intent(inout)  af,
real(dp), dimension(size(a, 1), size(a, 2)), intent(out)  b,
real(dp), dimension(size(af)), intent(out)  bf 
)
private

Definition at line 1647 of file mo_sce.F90.

◆ getpnt()

subroutine mo_sce::getpnt ( integer(i4), intent(in)  idist,
real(dp), dimension(:), intent(in)  bl,
real(dp), dimension(:), intent(in)  bu,
real(dp), dimension(:), intent(in)  std,
real(dp), dimension(:), intent(in)  xi,
logical, dimension(:), intent(in)  mask,
integer(i8), dimension(n_save_state), intent(inout)  save_state,
real(dp), dimension(size(xi, 1)), intent(out)  x 
)
private

Definition at line 1776 of file mo_sce.F90.

◆ parstt()

subroutine mo_sce::parstt ( real(dp), dimension(:, :), intent(in)  x,
real(dp), dimension(:), intent(in)  bound,
real(dp), intent(in)  peps,
logical, dimension(:), intent(in)  mask,
real(dp), dimension(size(bound, 1)), intent(out)  xnstd,
real(dp), intent(out)  gnrng,
integer(i4), intent(out)  ipcnvg 
)
private

Definition at line 1579 of file mo_sce.F90.

◆ sce()

real(dp) function, dimension(size(pini, 1)), public mo_sce::sce ( procedure(eval_interface), intent(in), pointer  eval,
procedure(objective_interface), intent(in), pointer  functn,
real(dp), dimension(:), intent(in)  pini,
real(dp), dimension(:, :), intent(in)  prange,
integer(i8), intent(in), optional  mymaxn,
logical, intent(in), optional  mymaxit,
integer(i4), intent(in), optional  mykstop,
real(dp), intent(in), optional  mypcento,
real(dp), intent(in), optional  mypeps,
integer(i8), intent(in), optional  myseed,
integer(i4), intent(in), optional  myngs,
integer(i4), intent(in), optional  mynpg,
integer(i4), intent(in), optional  mynps,
integer(i4), intent(in), optional  mynspl,
integer(i4), intent(in), optional  mymings,
integer(i4), intent(in), optional  myiniflg,
integer(i4), intent(in), optional  myprint,
logical, dimension(size(pini, 1)), intent(in), optional  mymask,
real(dp), intent(in), optional  myalpha,
real(dp), intent(in), optional  mybeta,
character(len = *), intent(in), optional  tmp_file,
character(len = *), intent(in), optional  popul_file,
logical, intent(in), optional  popul_file_append,
logical, intent(in), optional  parallel,
logical, intent(in), optional  restart,
character(len = *), intent(in), optional  restart_file,
real(dp), intent(out), optional  bestf,
integer(i8), intent(out), optional  neval,
real(dp), dimension(:), intent(out), optional, allocatable  history 
)

Shuffled Complex Evolution (SCE) algorithm for global optimization.

Shuffled Complex Evolution method for global optimization – version 2.1

by Qingyun Duan
Department of Hydrology & Water Resources
University of Arizona, Tucson, AZ 85721
(602) 621-9360, email: duan@.nosp@m.hwr..nosp@m.arizo.nosp@m.na.e.nosp@m.du

Written by Qingyun Duan, Oct 1990.
Revised by Qingyun Duan, Aug 1991.
Revised by Qingyun Duan, Apr 1992.

Re-written by Juliane Mai, Feb 2013.

Statement by Qingyun Duan:

This general purpose global optimization program is developed at the Department of Hydrology & Water Resources of the University of Arizona. Further information regarding the SCE-UA method can be obtained from Dr. Q. Duan, Dr. S. Sorooshian or Dr. V.K. Gupta at the address and phone number listed above. We request all users of this program make proper reference to the paper entitled 'Effective and Efficient Global Optimization for Conceptual Rainfall-runoff Models' by Duan, Q., S. Sorooshian, and V.K. Gupta, Water Resources Research, Vol 28(4), pp.1015-1031, 1992.

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

function functn(p)
use mo_kind, only: dp
implicit none
real(dp), dimension(:), intent(in) :: p
real(dp) :: functn
end function functn
Define number representations.
Definition mo_kind.F90:17
integer, parameter dp
Double Precision Real Kind.
Definition mo_kind.F90:46

Example

use mo_opt_functions, only: griewank
prange(:,1) = (/ -600.0, -600.0, -600.0, -600.0, -600.0, -600.0, -600.0, -600.0, -600.0, -600.0 /)
prange(:,2) = (/ 600.0, 600.0, 600.0, 600.0, 600.0, 600.0, 600.0, 600.0, 600.0, 600.0 /)
pini = (/ -.226265e+01, -.130187e+01, -.151219e+01, 0.133983e+00, 0.988159e+00, &
-.495074e+01, -.126574e+02, 0.572684e+00, 0.303864e+01, 0.343031e+01 /)
opt = sce(griewank, pini, prange)
Added for testing purposes of test_mo_sce, test_mo_dds, test_mo_mcmc.

-> see also example in test directory

Literature

  1. Duan, Q., S. Sorooshian, and V.K. Gupta - Effective and Efficient Global Optimization for Conceptual Rainfall-runoff Models Water Resources Research, Vol 28(4), pp.1015-1031, 1992.
Parameters
[in]real(dp) :: functn(p)Function on which to search the optimum
[in]real(dp) :: pini(:)inital value of decision variables
[in]real(dp) :: prange(size(pini),2)Min/max range of decision variables
[in]integer(i8), optional :: mymaxnmax no. of trials allowed before optimization is terminated
DEFAULT: 1000_i8
[in]logical, optional :: mymaxitmaximization (.true.) or minimization (.false.) of function
DEFAULT: false
[in]integer(i4), optional :: mykstopnumber of shuffling loops in which the criterion value must change by given percentage before optimiz. is terminated
DEFAULT: 10_i4
[in]real(dp), optional :: mypcentopercentage by which the criterion value must change in given number of shuffling loops
DEFAULT: 0.0001_dp
[in]real(dp), optional :: mypepsoptimization is terminated if volume of complex has converged to given percentage of feasible space
DEFAULT: 0.001_dp
[in]integer(i8), optional :: myseedinitial random seed
DEFAULT: get_timeseed
[in]integer(i4), optional :: myngsnumber of complexes in the initial population
DEFAULT: 2_i4
[in]integer(i4), optional :: mynpgnumber of points in each complex
DEFAULT: 2*n+1
[in]integer(i4), optional :: mynpsnumber of points in a sub-complex
DEFAULT: n+1
[in]integer(i4), optional :: mynsplnumber of evolution steps allowed for each complex before complex shuffling
DEFAULT: 2*n+1
[in]integer(i4), optional :: mymingsminimum number of complexes required, if the number of complexes is allowed to reduce as the optimization proceeds
DEFAULT: ngs = number of complexes in initial population
[in]integer(i4), optional :: myiniflgflag on whether to include the initial point in population
0, not included
1, included (DEFAULT)
[in]integer(i4), optional :: myprintflag for controlling print-out after each shuffling loop
0, print information on the best point of the population
1, print information on every point of the population
2, no printing (DEFAULT) 3, same as 0 but print progress '.' on every function call 4, same as 1 but print progress '.' on every function call
[in]logical, optional :: mymask(size(pini))parameter included in optimization (true) or discarded (false) DEFAULT: .true.
[in]real(dp), optional :: myalphaparameter for reflection of points in complex
DEFAULT: 0.8_dp
[in]real(dp), optional :: mybetaparameter for contraction of points in complex
DEFAULT: 0.45_dp
[in]character(len=*), optional :: tmp_fileif given: write results after each evolution loop to temporal output file of that name
# of headlines: 7
format: '# nloop icall ngs1 bestf worstf ...
... gnrng (bestx(j),j=1,nn)'
[in]character(len=*), optional :: popul_fileif given: write whole population to file of that name
# of headlines: 1
format: #_evolution_loop, xf(i), (x(i,j),j=1,nn)
total number of lines written <= neval <= mymaxn
[in]logical, optional :: popul_file_appendif true, append to existing population file (default: false)
[in]logical, optional :: parallelsce runs in parallel (true) or not (false) parallel sce should only be used if model/ objective is not parallel DEFAULT: .false.
[in]logical, optional :: restartif .true.: restart former sce run from restart_file
[in]character(len=*), optional :: restart_filefile name for read/write of restart file (default: mo_sce.restart)
[out]real(dp), optional :: bestfthe best value of the function.
[out]integer(i8), optional :: nevalnumber of function evaluations needed.
[out]real(dp), optional, allocatable :: history(:)the history of best function values, history(neval)=bestf
Returns
real(dp) :: bestx(size(pini)) — The parameters of the point which is estimated to minimize/maximize the function.
Note
Maximal number of parameters is 1000.
SCE is OpenMP enabled on the loop over the complexes.
OMP_NUM_THREADS > 1 does not give reproducible results even when seeded!
Author
Juliane Mai, Matthias Cuntz
Date
Feb 2013
  • first version
Jul 2013
  • OpenMP
  • NaN and Inf in objective function
Oct 2013
  • added peps as optional argument
  • allow for masked parameters
  • write population to file --> popul_file
  • write intermediate results to file --> tmp_file
  • flag parallel introduced
Nov 2013
  • progress dots (MC)
  • use iso_fortran_env
  • treat functn=NaN as worse function value in cce
May 2014
  • sort -> orderpack (MC)
  • popul_file_append (MC)
  • sort with NaNs (MC)
Feb 2015
  • restart (MC, JM)
Mar 2015
  • use is_finite and special_value from mo_utils (MC)
Apr 2015
  • handling of array-like variables in restart-namelists (JM)

Definition at line 212 of file mo_sce.F90.

References mo_string_utils::compress(), mo_kind::dp, mo_kind::i4, mo_kind::i8, and mo_xor4096::n_save_state.

Here is the call graph for this function:

◆ sort_matrix()

subroutine mo_sce::sort_matrix ( real(dp), dimension(:, :), intent(inout)  rb,
real(dp), dimension(:), intent(inout)  ra 
)
private

Definition at line 1703 of file mo_sce.F90.