initialize_sa Subroutine

private subroutine initialize_sa(me, fcn, n, lb, ub, c, maximize, eps, ns, nt, neps, maxevl, iprint, iseed1, iseed2, step_mode, vms, iunit, use_initial_guess, n_resets, optimal_f_specified, optimal_f, optimal_f_tol)

Initialize the class.

See the definition of simulated_annealing_type for the options. Any options not set here will use the default values in the type.

Type Bound

simulated_annealing_type

Arguments

Type IntentOptional Attributes Name
class(simulated_annealing_type), intent(inout) :: me
procedure(sa_func) :: fcn
integer, intent(in) :: n
real(kind=wp), intent(in), dimension(n) :: lb
real(kind=wp), intent(in), dimension(n) :: ub
real(kind=wp), intent(in), optional, dimension(n) :: c
logical, intent(in), optional :: maximize
real(kind=wp), intent(in), optional :: eps
integer, intent(in), optional :: ns
integer, intent(in), optional :: nt
integer, intent(in), optional :: neps
integer, intent(in), optional :: maxevl
integer, intent(in), optional :: iprint
integer, intent(in), optional :: iseed1
integer, intent(in), optional :: iseed2
integer, intent(in), optional :: step_mode
real(kind=wp), intent(in), optional :: vms
integer, intent(in), optional :: iunit
logical, intent(in), optional :: use_initial_guess
integer, intent(in), optional :: n_resets
logical, intent(in), optional :: optimal_f_specified

if the optional f value is known, it can be specified by optimal_f. [Default is False]

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

if optimal_f_specified=True the solver will stop if this value is achieved.

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

absolute tolerance for the optimal_f check


Calls

proc~~initialize_sa~~CallsGraph proc~initialize_sa simulated_annealing_module::simulated_annealing_type%initialize_sa proc~destroy_sa simulated_annealing_module::simulated_annealing_type%destroy_sa proc~initialize_sa->proc~destroy_sa

Source Code

    subroutine initialize_sa(me,fcn,n,lb,ub,c,&
                             maximize,eps,ns,nt,neps,maxevl,&
                             iprint,iseed1,iseed2,step_mode,vms,iunit,&
                             use_initial_guess,n_resets,&
                             optimal_f_specified,optimal_f,optimal_f_tol)

    implicit none

    class(simulated_annealing_type),intent(inout) :: me
    procedure(sa_func)                            :: fcn
    integer, intent(in)                           :: n
    real(wp), dimension(n), intent(in)            :: lb
    real(wp), dimension(n), intent(in)            :: ub
    real(wp), dimension(n), intent(in),optional   :: c
    logical, intent(in),optional                  :: maximize
    real(wp), intent(in),optional                 :: eps
    integer, intent(in),optional                  :: ns
    integer, intent(in),optional                  :: nt
    integer, intent(in),optional                  :: neps
    integer, intent(in),optional                  :: maxevl
    integer, intent(in),optional                  :: iprint
    integer, intent(in),optional                  :: iseed1
    integer, intent(in),optional                  :: iseed2
    integer,intent(in),optional                   :: step_mode
    real(wp), intent(in), optional                :: vms
    integer, intent(in), optional                 :: iunit
    logical, intent(in), optional                 :: use_initial_guess
    integer, intent(in), optional                 :: n_resets
    logical, intent(in), optional  :: optimal_f_specified !! if the optional `f` value is known,
                                                          !! it can be specified by `optimal_f`.
                                                          !! [Default is False]
    real(wp), intent(in), optional  :: optimal_f          !! if `optimal_f_specified=True` the solver will stop
                                                          !! if this value is achieved.
    real(wp), intent(in), optional  :: optimal_f_tol      !! absolute tolerance for the `optimal_f` check

    call me%destroy()

    me%fcn => fcn
    me%n  = n
    me%lb = lb
    me%ub = ub

    ! check validity of bounds:
    if (any(me%lb > me%ub)) then
        error stop 'Error: at least one LB > UB.'
    end if

    allocate(me%c(n))
    if (present(c)) then
        me%c = c
    else
        me%c = 2.0_wp
    end if

    if (present(maximize)            ) me%maximize = maximize
    if (present(eps)                 ) me%eps = abs(eps)
    if (present(ns)                  ) me%ns = ns
    if (present(nt)                  ) me%nt = nt
    if (present(neps)                ) me%neps = neps
    if (present(maxevl)              ) me%maxevl = maxevl
    if (present(iprint)              ) me%iprint = iprint
    if (present(iseed1)              ) me%iseed1 = iseed1
    if (present(iseed2)              ) me%iseed2 = iseed2
    if (present(step_mode)           ) me%step_mode = step_mode
    if (present(vms)                 ) me%vms = vms
    if (present(iunit)               ) me%iunit = iunit
    if (present(use_initial_guess)   ) me%use_initial_guess = use_initial_guess
    if (present(n_resets)            ) me%n_resets = n_resets
    if (present(optimal_f_specified) ) me%optimal_f_specified = optimal_f_specified
    if (present(optimal_f)           ) me%optimal_f = optimal_f
    if (present(optimal_f_tol)       ) me%optimal_f_tol = abs(optimal_f_tol)

    end subroutine initialize_sa