set_inputs Subroutine

private subroutine set_inputs(me, n, xl, xu, f, status, iter_f, np, ngen, nd, pcross, pmutmn, pmutmx, pmut, imut, fdif, irep, ielite, ivrb, convergence_tol, convergence_window, initial_guess_frac, iseed)

Constructor for the pikaia_class. The routine must be called before the solve routine can be used.

The following inputs are required: n, f, xl, xu. For the others, if they are not present, then the default values are used

Note

Based on setctl in the original code.

Type Bound

pikaia_class

Arguments

Type IntentOptional Attributes Name
class(pikaia_class), intent(out) :: me

pikaia class

integer, intent(in) :: n

the parameter space dimension, i.e., the number of adjustable parameters (size of the x vector).

real(kind=wp), intent(in), dimension(n) :: xl

vector of lower bounds for x

real(kind=wp), intent(in), dimension(n) :: xu

vector of upper bounds for x

procedure(pikaia_func) :: f

user-supplied scalar function of n variables, which must have the pikaia_func procedure interface. By convention, f should return higher values for more optimal parameter values (i.e., individuals which are more "fit"). For example, in fitting a function through data points, f could return the inverse of chi**2.

integer, intent(out) :: status

status output flag (0 if there were no errors)

procedure(iter_func), optional :: iter_f

user-supplied subroutine that will report the best solution for each generation. It must have the iter_func procedure interface. If not present, then it is not used. (note: this is independent of ivrb).

integer, intent(in), optional :: np

number of individuals in a population (default is 100)

integer, intent(in), optional :: ngen

maximum number of iterations

integer, intent(in), optional :: nd

number of significant digits (i.e., number of genes) retained in chromosomal encoding (default is 6).

real(kind=wp), intent(in), optional :: pcross

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

real(kind=wp), intent(in), optional :: pmutmn

minimum mutation rate; must be >= 0.0 (default is 0.0005)

real(kind=wp), intent(in), optional :: pmutmx

maximum mutation rate; must be <= 1.0 (default is 0.25)

real(kind=wp), intent(in), optional :: pmut

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.)

integer, intent(in), optional :: imut

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.

real(kind=wp), intent(in), optional :: fdif

relative fitness differential; range from 0 (none) to 1 (maximum). (default is 1.0)

integer, intent(in), optional :: irep

reproduction plan; 1/2/3=Full generational replacement/Steady-state-replace-random/Steady- state-replace-worst (default is 3)

integer, intent(in), optional :: ielite

elitism flag; 0/1=off/on (default is 0) (Applies only to reproduction plans 1 and 2)

integer, intent(in), optional :: ivrb

printed output 0/1/2=None/Minimal/Verbose (default is 0)

real(kind=wp), intent(in), optional :: convergence_tol

convergence tolerance; must be > 0.0 (default is 0.0001)

integer, intent(in), optional :: convergence_window

convergence window; must be >= 0 This is the number of consecutive solutions within the tolerance for convergence to be declared (default is 20)

real(kind=wp), intent(in), optional :: initial_guess_frac

fraction 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%).

integer, intent(in), optional :: iseed

random seed value; must be > 0 (default is 999)


Source Code

    subroutine set_inputs(me,&
                            n,xl,xu,f,status,&
                            iter_f,&
                            np,ngen,nd,pcross,pmutmn,pmutmx,pmut,imut,&
                            fdif,irep,ielite,ivrb,&
                            convergence_tol,convergence_window,initial_guess_frac,&
                            iseed)

    implicit none

    class(pikaia_class),intent(out)    :: me        !! pikaia class
    integer,intent(in)                 :: n         !! the parameter space dimension, i.e., the number
                                                    !! of adjustable parameters (size of the x vector).
    real(wp),dimension(n),intent(in)   :: xl        !! vector of lower bounds for x
    real(wp),dimension(n),intent(in)   :: xu        !! vector of upper bounds for x
    procedure(pikaia_func)             :: f         !! user-supplied scalar function of n variables,
                                                    !! which must have the [[pikaia_func]] procedure interface.
                                                    !! By convention, f should return higher values for more optimal
                                                    !! parameter values (i.e., individuals which are more "fit").
                                                    !! For example, in fitting a function through data points, f
                                                    !! could return the inverse of chi**2.
    integer,intent(out)                :: status    !! status output flag (0 if there were no errors)
    procedure(iter_func),optional      :: iter_f    !! user-supplied subroutine that will report the
                                                    !! best solution for each generation.
                                                    !! It must have the [[iter_func]] procedure interface.  If not present,
                                                    !! then it is not used.  (note: this is independent of ivrb).
    integer,intent(in),optional        :: np        !! number of individuals in a population (default is 100)
    integer,intent(in),optional        :: ngen      !! maximum number of iterations
    integer,intent(in),optional        :: nd        !! number of significant digits (i.e., number of
                                                    !! genes) retained in chromosomal encoding (default is 6).
    real(wp),intent(in),optional       :: pcross    !! 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
    real(wp),intent(in),optional       :: pmutmn    !! minimum mutation rate; must be >= 0.0 (default is 0.0005)
    real(wp),intent(in),optional       :: pmutmx    !! maximum mutation rate; must be <= 1.0 (default is 0.25)
    real(wp),intent(in),optional       :: pmut      !! 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.)
    integer,intent(in),optional        :: imut      !! 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.
    real(wp),intent(in),optional       :: fdif      !! relative fitness differential; range from 0
                                                    !! (none) to 1 (maximum).  (default is 1.0)
    integer,intent(in),optional        :: irep      !! reproduction plan; 1/2/3=Full generational
                                                    !! replacement/Steady-state-replace-random/Steady-
                                                    !! state-replace-worst (default is 3)
    integer,intent(in),optional        :: ielite    !! elitism flag; 0/1=off/on (default is 0)
                                                    !! (Applies only to reproduction plans 1 and 2)
    integer,intent(in),optional        :: ivrb      !! printed output 0/1/2=None/Minimal/Verbose
                                                    !! (default is 0)
    real(wp),intent(in),optional       :: convergence_tol    !! convergence tolerance; must be > 0.0 (default is 0.0001)
    integer,intent(in),optional        :: convergence_window !! convergence window; must be >= 0
                                                             !! This is the number of consecutive solutions
                                                             !! within the tolerance for convergence to
                                                             !! be declared (default is 20)
    real(wp),intent(in),optional       :: initial_guess_frac !! fraction 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%).
    integer,intent(in),optional        :: iseed              !! random seed value; must be > 0 (default is 999)

    me%n = n

    if (allocated(me%xl)) deallocate(me%xl)
    allocate(me%xl(n))
    me%xl = xl

    if (allocated(me%xu)) deallocate(me%xu)
    allocate(me%xu(n))
    me%xu = xu

    if (allocated(me%del)) deallocate(me%del)
    allocate(me%del(n))
    me%del = me%xu - me%xl

    me%user_f => f

    if (present(iter_f)) me%iter_f => iter_f

    if (present(np                 )) me%np                 = np
    if (present(ngen               )) me%ngen               = ngen
    if (present(nd                 )) me%nd                 = nd
    if (present(pcross             )) me%pcross             = pcross
    if (present(imut               )) me%imut               = imut
    if (present(pmut               )) me%pmuti              = pmut  !initial value
    if (present(pmutmn             )) me%pmutmn             = pmutmn
    if (present(pmutmx             )) me%pmutmx             = pmutmx
    if (present(fdif               )) me%fdif               = fdif
    if (present(irep               )) me%irep               = irep
    if (present(ielite             )) me%ielite             = ielite
    if (present(ivrb               )) me%ivrb               = ivrb
    if (present(convergence_tol    )) me%convergence_tol    = convergence_tol
    if (present(convergence_window )) me%convergence_window = convergence_window
    if (present(initial_guess_frac )) me%initial_guess_frac = initial_guess_frac
    if (present(iseed              )) me%iseed              = iseed

    !check for errors:

    !initialize error flag:
    status = 0

    !Print a header
    if (me%ivrb>0) then
        write(output_unit,'(A)') '------------------------------------------------------------'
        write(output_unit,'(A)') '              PIKAIA Genetic Algorithm Report               '
        write(output_unit,'(A)') '------------------------------------------------------------'
        write(output_unit,'(A,I4)')    ' Number of Generations evolving: ',me%ngen
        write(output_unit,'(A,I4)')    '     Individuals per generation: ',me%np
        write(output_unit,'(A,I4)')    '  Number of Chromosome segments: ',me%n
        write(output_unit,'(A,I4)')    '  Length of Chromosome segments: ',me%nd
        write(output_unit,'(A,E11.4)') '          Crossover probability: ',me%pcross
        write(output_unit,'(A,E11.4)') '          Initial mutation rate: ',me%pmuti
        write(output_unit,'(A,E11.4)') '          Minimum mutation rate: ',me%pmutmn
        write(output_unit,'(A,E11.4)') '          Maximum mutation rate: ',me%pmutmx
        write(output_unit,'(A,E11.4)') '  Relative fitness differential: ',me%fdif
        write(output_unit,'(A,E11.4)') '         Initial guess fraction: ',me%initial_guess_frac
        write(output_unit,'(A,E11.4)') '          Convergence tolerance: ',me%convergence_tol
        write(output_unit,'(A,I4)')    '             Convergence window: ',me%convergence_window
        select case (me%imut)
        case(1); write(output_unit,'(A)') '                  Mutation Mode: Uniform, Constant Rate'
        case(2); write(output_unit,'(A)') '                  Mutation Mode: Uniform, Variable Rate (F)'
        case(3); write(output_unit,'(A)') '                  Mutation Mode: Uniform, Variable Rate (D)'
        case(4); write(output_unit,'(A)') '                  Mutation Mode: Uniform+Creep, Constant Rate'
        case(5); write(output_unit,'(A)') '                  Mutation Mode: Uniform+Creep, Variable Rate (F)'
        case(6); write(output_unit,'(A)') '                  Mutation Mode: Uniform+Creep, Variable Rate (D)'
        end select
        select case (me%irep)
        case(1); write(output_unit,'(A)') '              Reproduction Plan: Full generational replacement'
        case(2); write(output_unit,'(A)') '              Reproduction Plan: Steady-state-replace-random'
        case(3); write(output_unit,'(A)') '              Reproduction Plan: Steady-state-replace-worst'
        end select
        write(output_unit,'(A)') '------------------------------------------------------------'
    end if

    !Check some control values
    if (me%imut/=1 .and. me%imut/=2 .and. me%imut/=3 .and. &
          me%imut/=4 .and. me%imut/=5 .and. me%imut/=6) then
       write(output_unit,'(A)') ' ERROR: illegal value for Mutation Mode.'
       status = 5
    end if

    if (me%fdif>1.0_wp) then
       write(output_unit,'(A)') ' ERROR: illegal value for Relative fitness differential.'
       status = 9
    end if

    if (me%irep/=1 .and. me%irep/=2 .and. me%irep/=3) then
       write(output_unit,'(A)') ' ERROR: illegal value for Reproduction plan.'
       status = 10
    end if

    if (me%pcross>1.0_wp .or. me%pcross<0.0_wp) then
       write(output_unit,'(A)') ' ERROR: illegal value for Crossover probability.'
       status = 4
    end if

    if (me%ielite/=0 .and. me%ielite/=1) then
       write(output_unit,'(A)') ' ERROR: illegal value for Elitism flag.'
       status = 11
    end if

    if (me%convergence_tol<=0.0_wp) then
        write(output_unit,'(A)') ' ERROR: illegal value for Convergence tolerance.'
        status = 101
    end if

    if (me%convergence_window<=0) then
        write(output_unit,'(A)') ' ERROR: illegal value for Convergence window.'
        status = 102
    end if

    if (me%iseed<=0) then
        write(output_unit,'(A)') ' ERROR: illegal value for iseed.'
        status = 103
    end if

    if (me%nd>9 .or. me%nd<1) then
        write(output_unit,'(A)') ' ERROR: illegal value for Chromosome length.'
        status = 104
    end if

    if (mod(me%np,2)>0) then
       write(output_unit,'(A)') ' ERROR: population size must be an even number.'
       status = 105
    end if

    if (me%initial_guess_frac<0.0_wp .or. me%initial_guess_frac>1.0_wp) then
       write(output_unit,'(A)') ' ERROR: illegal value for Initial guess fraction.'
       status = 106
    end if

    if (me%irep==1 .and. me%imut==1 .and. me%pmuti>0.5_wp .and. me%ielite==0) then
       write(output_unit,'(A)') &
        ' WARNING: dangerously high value for Initial mutation rate; '//&
        '(Should enforce elitism with ielite=1.)'
    end if

    if (me%irep==1 .and. me%imut==2 .and. me%pmutmx>0.5_wp .and. me%ielite==0) then
       write(output_unit,'(A)') &
       ' WARNING: dangerously high value for Maximum mutation rate; '//&
       '(Should enforce elitism with ielite=1.)'
    end if

    if (me%fdif<0.33_wp .and. me%irep/=3) then
       write(output_unit,'(A)') &
       ' WARNING: dangerously low value of Relative fitness differential.'
    end if

    end subroutine set_inputs