Note
This is a refactoring of the PIKAIA unconstrained optimization code from the High Altitude Observatory. The original code is public domain and was written by Paul Charbonneau & Barry Knapp.
The present code is the awesome modern Fortran version written by Jabob Williams:
Type | Intent | Optional | Attributes | Name | ||||||||||||||||||||||||||||||
---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|
type(pikaia_class), | intent(inout) | :: | pik_class |
PIKAIA class instanciation |
||||||||||||||||||||||||||||||
character(len=4), | intent(in) | :: | step |
init or solv |
||||||||||||||||||||||||||||||
real(kind=R8), | intent(in), | dimension(1:nparam) | :: | xl |
lower bonds of xx |
|||||||||||||||||||||||||||||
real(kind=R8), | intent(in), | dimension(1:nparam) | :: | xu |
upper bonds of xx |
|||||||||||||||||||||||||||||
integer(kind=I4), | intent(in) | :: | nparam |
number of parameters |
||||||||||||||||||||||||||||||
private subroutine cost(me, x, f)Arguments
|
||||||||||||||||||||||||||||||||||
integer(kind=I4), | intent(out) | :: | istat | |||||||||||||||||||||||||||||||
real(kind=R8), | intent(out) | :: | f | |||||||||||||||||||||||||||||||
real(kind=R8), | intent(out), | dimension(1:nparam) | :: | xx |
chromosom for PIKAIA |
subroutine pikaia_skku_solver(pik_class, step, xl, xu, nparam, cost, istat, f, xx) !================================================================================================ !<@note This is a refactoring of the PIKAIA unconstrained optimization code from the High Altitude Observatory. !< The original code is public domain and was written by Paul Charbonneau & Barry Knapp. !< !< The present code is the awesome modern Fortran version written by Jabob Williams: !< !< [OOP Pikaia, Jacob Williams](https://github.com/jacobwilliams/pikaia) !< !<@endnote !------------------------------------------------------------------------------------------------ implicit none type(pikaia_class), intent(inout) :: pik_class !! **PIKAIA** *class instanciation* character(len=4), intent(in ) :: step !! *init* or *solv* integer(kind=I4), intent(in ) :: nparam !! *number of parameters* real(kind=R8), intent(in ), dimension(1:nparam) :: xl !! *lower bonds of xx* real(kind=R8), intent(in ), dimension(1:nparam) :: xu !! *upper bonds of xx* real(kind=R8), intent( out), dimension(1:nparam) :: xx !! *chromosom for* **PIKAIA** integer(kind=I4), intent( out) :: istat real(kind=R8), intent( out) :: f interface subroutine cost(me, x, f) use data_arch, only : R8 use pikaia_oop, only : pikaia_class implicit none class(pikaia_class), intent(inout) :: me real(kind=R8) , intent(in ), dimension(:) :: x real(kind=R8) , intent( out) :: f endsubroutine cost endinterface select case(step) case('init') !initialize the class: call pik_class%init( n = nparam, & ! IN ; the parameter space dimension, i.e., the number of adjustable parameters (size of the x vector). xl = xl, & ! IN, DIM(n) ; vector of lower bounds for x xu = xu, & ! IN, DIM(n) ; vector of upper bounds for x f = cost, & ! ; user-supplied scalar function of n variables, which must have the pikaia_func procedure interface. status = istat , & ! OUT ; status output flag (0 if there were no errors) !iter_f = report_iteration, & ! OPT ; user-supplied subroutine that will report the best solution for each generation. It must have the iter_func procedure interface. np = 100, & ! IN, OPT ; number of individuals in a population (default is 100) ngen = 1000, & ! IN, OPT ; maximum number of iterations nd = 9, & ! IN ; number of significant digits (i.e., number of genes) retained in chromosomal encoding pcross = 0.85_R8, & ! IN, OPT ; crossover probability; must be <= 1.0 (default is 0.85). If crossover takes place, either one or two splicing points are used, with equal probabilities pmutmn = 0.0005_R8, & ! IN, OPT ; minimum mutation rate; must be >= 0.0 (default is 0.0005) pmutmx = 0.25_R8, & ! IN, OPT ; maximum mutation rate; must be <= 1.0 (default is 0.25) pmut = 0.005_R8, & ! IN, OPT ; initial mutation rate; should be small (default is 0.005) (Note: the mutation rate is the probability that any one gene locus will mutate in any one generation.) imut = 2, & ! IN, OPT ; mutation mode; 1/2/3/4/5 (default is 2). ! 1=one-point mutation, fixed rate. ! 2=one-point, adjustable rate based on fitness. ! 3=one-point, adjustable rate based on distance. ! 4=one-point+creep, fixed rate. ! 5=one-point+creep, adjustable rate based on fitness. ! 6=one-point+creep, adjustable rate based on distance. fdif = 1._R8, & ! IN, OPT ; relative fitness differential; range from 0 (none) to 1 (maximum). (default is 1.0) irep = 3, & ! IN, OPT ; reproduction plan; 1/2/3=Full generational replacement/Steady-state-replace-random/Steady- state-replace-worst (default is 3) ielite = 0, & ! IN, OPT ; elitism flag; 0/1=off/on (default is 0) (Applies only to reproduction plans 1 and 2) ivrb = 0, & ! IN, OPT ; printed output 0/1/2=None/Minimal/Verbose convergence_tol = 1.0e-6_R8, & ! IN, OPT ; convergence tolerance; must be > 0.0 (default is 0.0001) convergence_window = 200, & ! IN, OPT ; convergence window; must be >= 0 This is the number of consecutive solutions within the tolerance for convergence to be declared (default is 20) initial_guess_frac = 0.1_R8, & ! IN, OPT ; raction of the initial population to set equal to the initial guess. Range from 0 (none) to 1.0 (all). (default is 0.1 or 10%). iseed = 999) ! IN, OPT ; random seed value; must be > 0 (default is 999) case('solv') call pik_class%solve( x = xx(1:nparam), & ! INOUT, DIM(*) ; f = f, & ! OUT ; status = istat, & ! OUT ; omp = .true. ) ! IN, OPTIONAL case default stop 'Wrong choice in "pikaia_skku_solver"' endselect return endsubroutine pikaia_skku_solver