slsqp_wrapper Subroutine

private subroutine slsqp_wrapper(me, x, istat, iterations, status_message)

main routine for calling slsqp.

Type Bound

slsqp_solver

Arguments

Type IntentOptional Attributes Name
class(slsqp_solver), intent(inout) :: me
real(kind=wp), intent(inout), dimension(:) :: x

in: initial optimization variables, out: solution.

integer, intent(out) :: istat

status code (see mode in slsqp).

integer, intent(out), optional :: iterations

number of iterations

character(len=:), intent(out), optional, allocatable :: status_message

string status message corresponding to istat


Calls

proc~~slsqp_wrapper~~CallsGraph proc~slsqp_wrapper slsqp_module::slsqp_solver%slsqp_wrapper proc~destroy_linmin_data slsqp_core::linmin_data%destroy_linmin_data proc~slsqp_wrapper->proc~destroy_linmin_data proc~mode_to_status_message slsqp_module::mode_to_status_message proc~slsqp_wrapper->proc~mode_to_status_message proc~report_message slsqp_module::slsqp_solver%report_message proc~slsqp_wrapper->proc~report_message proc~slsqp slsqp_core::slsqp proc~slsqp_wrapper->proc~slsqp proc~slsqpb slsqp_core::slsqpb proc~slsqp->proc~slsqpb proc~check_convergence slsqp_core::check_convergence proc~slsqpb->proc~check_convergence proc~daxpy slsqp_support::daxpy proc~slsqpb->proc~daxpy proc~dcopy slsqp_support::dcopy proc~slsqpb->proc~dcopy proc~ddot slsqp_support::ddot proc~slsqpb->proc~ddot proc~dscal slsqp_support::dscal proc~slsqpb->proc~dscal proc~enforce_bounds slsqp_core::enforce_bounds proc~slsqpb->proc~enforce_bounds proc~ldl slsqp_core::ldl proc~slsqpb->proc~ldl proc~linmin slsqp_core::linmin proc~slsqpb->proc~linmin proc~lsq slsqp_core::lsq proc~slsqpb->proc~lsq proc~dnrm2 slsqp_support::dnrm2 proc~check_convergence->proc~dnrm2 proc~lsq->proc~dcopy proc~lsq->proc~ddot proc~lsq->proc~dscal proc~lsq->proc~enforce_bounds proc~lsei slsqp_core::lsei proc~lsq->proc~lsei proc~lsei->proc~dcopy proc~lsei->proc~ddot proc~lsei->proc~dnrm2 proc~h12 slsqp_core::h12 proc~lsei->proc~h12 proc~hfti slsqp_core::hfti proc~lsei->proc~hfti proc~lsi slsqp_core::lsi proc~lsei->proc~lsi proc~hfti->proc~h12 proc~lsi->proc~daxpy proc~lsi->proc~ddot proc~lsi->proc~dnrm2 proc~lsi->proc~h12 proc~ldp slsqp_core::ldp proc~lsi->proc~ldp proc~ldp->proc~daxpy proc~ldp->proc~dcopy proc~ldp->proc~ddot proc~ldp->proc~dnrm2 proc~bvls_wrapper bvls_module::bvls_wrapper proc~ldp->proc~bvls_wrapper proc~nnls slsqp_core::nnls proc~ldp->proc~nnls proc~bvls bvls_module::bvls proc~bvls_wrapper->proc~bvls proc~nnls->proc~h12 proc~g1 slsqp_core::g1 proc~nnls->proc~g1

Source Code

    subroutine slsqp_wrapper(me,x,istat,iterations,status_message)

    implicit none

    class(slsqp_solver),intent(inout)   :: me
    real(wp),dimension(:),intent(inout) :: x          !! **in:**  initial optimization variables,
                                                      !! **out:** solution.
    integer,intent(out)                 :: istat      !! status code (see `mode` in [[slsqp]]).
    integer,intent(out),optional        :: iterations !! number of iterations
    character(len=:),intent(out),allocatable,optional :: status_message
                                                      !! string status message
                                                      !! corresponding to `istat`

    ! local variables:
    real(wp),dimension(:),allocatable   :: c        !! constraint vector -- `dimension(max(1,me%m))`
    real(wp),dimension(:,:),allocatable :: a        !! a matrix for [[slsqp]] -- `dimension(max(1,me%m),me%n+1)`
    real(wp),dimension(:),allocatable   :: g        !! g matrix for [[slsqp]] -- `dimension(me%n+1)`
    real(wp),dimension(:),allocatable   :: cvec     !! constraint vector -- `dimension(me%m)`
    real(wp),dimension(:),allocatable   :: dfdx     !! objective function partials -- `dimension(me%n)`
    real(wp),dimension(:,:),allocatable :: dcdx     !! constraint partials -- `dimension(me%m,me%n)`
    real(wp),dimension(:),allocatable   :: delta    !! perturbation step size to approximate gradient -- `dimension(me%n)`
    real(wp),dimension(:),allocatable   :: cvecr    !! right function value to approximate constraints vector's gradient -- `dimension(me%m)`
    real(wp),dimension(:),allocatable   :: cvecl    !! left function value to approximate constraints vector's gradient -- `dimension(me%m)`
    real(wp)                            :: f        !! objective function
    integer                             :: i        !! iteration counter
    integer                             :: mode     !! reverse communication flag for [[slsqp]]
    integer                             :: la       !! input to [[slsqp]]
    integer                             :: iter     !! in/out for [[slsqp]]
    real(wp)                            :: acc      !! in/out for [[slsqp]]
    integer                             :: ig       !! loop index to approximate gradient
    real(wp)                            :: fr       !! right function value to approximate objective function's gradient
    real(wp)                            :: fl       !! left function value to approximate objective function's gradient
    real(wp)                            :: fact     !! denominator factor for finite difference approximation

    !initialize:
    allocate(c(max(1,me%m))       )
    allocate(a(max(1,me%m),me%n+1))
    allocate(g(me%n+1)            )
    allocate(cvec(me%m)           )
    allocate(dfdx(me%n)           )
    allocate(dcdx(me%m,me%n)      )
    allocate(delta(me%n)          )
    allocate(cvecr(me%m)          )
    allocate(cvecl(me%m)          )
    i    = 0
    iter = me%max_iter
    la   = max(1,me%m)
    mode = 0
    a    = zero
    g    = zero
    c    = zero
    if (present(iterations)) iterations = 0
    call me%linmin%destroy()
    call me%slsqpb%destroy()

    !check setup:
    if (size(x)/=me%n) then
        istat = -100
        call me%report_message(mode_to_status_message(istat))
        if (present(status_message)) status_message = mode_to_status_message(istat)
        return
    end if

    !linesearch:
    select case(me%linesearch_mode)
    case(1)     !inexact (armijo-type linesearch)
        acc = abs(me%acc)
    case(2)     !exact
        acc = -abs(me%acc)
    case default
        istat = -101
        call me%report_message(mode_to_status_message(istat))
        if (present(status_message)) status_message = mode_to_status_message(istat)
        return
    end select

    !make sure the functions have been associated:
    if (.not. associated(me%f)) then
        istat = -102
        call me%report_message(mode_to_status_message(istat))
        if (present(status_message)) status_message = mode_to_status_message(istat)
        return
    end if
    if ((me%gradient_mode==0).and.(.not. associated(me%g))) then
        istat = -103
        call me%report_message(mode_to_status_message(istat))
        if (present(status_message)) status_message = mode_to_status_message(istat)
        return
    end if
    if (me%gradient_mode<0 .or. me%gradient_mode>3) then
        istat = -104
        call me%report_message(mode_to_status_message(istat))
        if (present(status_message)) status_message = mode_to_status_message(istat)
        return
    end if
    if (me%gradient_mode/=0 .and. me%gradient_delta<=epmach) then
        istat = -105
        call me%report_message(mode_to_status_message(istat))
        if (present(status_message)) status_message = mode_to_status_message(istat)
        return
    end if

    !main solver loop:
    do

        if (mode==0 .or. mode==1) then  !function evaluation (f&c)
            call me%f(x,f,cvec)
            c(1:me%m) = cvec
        end if

        if (mode==0 .or. mode==-1) then  !gradient evaluation (g&a)
            select case (me%gradient_mode)
            case (0) ! user supplied gradients
                call me%g(x,dfdx,dcdx)
                g(1:me%n)        = dfdx
                a(1:me%m,1:me%n) = dcdx
            case default ! approximate using finite differences
                if (me%gradient_mode==3) then
                    fact = two  ! central differences
                else
                    fact = one  ! forward/backward differences
                end if
                do ig=1,me%n
                    !initialize a delta to perturb the objective
                    !function and the constraint vector
                    delta     = zero
                    delta(ig) = me%gradient_delta
                    !get the right and left value of the objective
                    !function and the constraint vector
                    select case (me%gradient_mode)
                    case (1) ! backward difference
                        call me%f(x,fr,cvecr)
                        call me%f(x-delta,fl,cvecl)
                    case (2) ! forward difference
                        call me%f(x+delta,fr,cvecr)
                        call me%f(x,fl,cvecl)
                    case (3) ! central difference
                        call me%f(x+delta,fr,cvecr)
                        call me%f(x-delta,fl,cvecl)
                    end select
                    !compute the gradients by first-order finite differences
                    g(ig) = (fr-fl) / ( fact*delta(ig) )
                    if (me%m>0) then
                        a(:,ig) = (cvecr-cvecl) / ( fact*delta(ig) )
                    end if
                end do
            end select
            !this is an iteration:
            !note: the initial guess is reported as iteration 0:
            if (associated(me%report)) call me%report(i,x,f,c) !report iteration
            i = i + 1  ! iteration counter
        end if

        !main routine:
        call slsqp(me%m,me%meq,la,me%n,x,me%xl,me%xu,&
                   f,c,g,a,acc,iter,mode,&
                   me%w,me%l_w, &
                   me%slsqpb,me%linmin,me%alphamin,me%alphamax,&
                   me%tolf,me%toldf,me%toldx,&
                   me%max_iter_ls,me%nnls_mode,&
                   me%infinite_bound)

        if (mode==1 .or. mode==-1) then
            !continue to next call
        else
            if (mode==0 .and. associated(me%report)) &
                call me%report(i,x,f,c) !report solution
            call me%report_message(mode_to_status_message(mode))
            exit
        end if

        if (me%user_triggered_stop) then
            mode = -2
            call me%report_message(mode_to_status_message(mode))
            me%user_triggered_stop = .false.    !have to reset in case
                                                !method is called again.
            exit
        end if

    end do

    istat = mode
    if (present(iterations)) iterations = iter
    if (present(status_message)) status_message = mode_to_status_message(istat)

    end subroutine slsqp_wrapper