psqpn Subroutine

private subroutine psqpn(me, nf, nb, nc, x, bound_type, xl, xu, cf, constraint_type, cl, cu, ipar, rpar, f, gmax, cmax, iprnt, iterm, obj, dobj, con, dcon, report)

easy to use subroutine for general nonlinear programming problems.

Type Bound

psqp_class

Arguments

Type IntentOptional Attributes Name
class(psqp_class), intent(inout) :: me
integer, intent(in) :: nf

number of variables

integer, intent(in) :: nb

choice of simple bounds.

  • nb=0 -- simple bounds suppressed.
  • nb>0 -- simple bounds accepted.
integer, intent(in) :: nc

number of general nonlinear constraints.

real(kind=wp), intent(inout), dimension(nf) :: x

x(nf) vector of variables.

integer, intent(in), dimension(nf) :: bound_type

ix(nf) vector containing types of bounds.

  • ix(i) = 0 -- variable x(i) is unbounded.
  • ix(i) = 1 -- lower bound xl(i) <= x(i).
  • ix(i) = 2 -- upper bound x(i) <= xu(i).
  • ix(i) = 3 -- two side bound xl(i) <= x(i) <= xu(i).
  • ix(i) = 5 -- variable x(i) is fixed.
real(kind=wp), intent(in), dimension(nf) :: xl

xl(nf) vector containing lower bounds for variables.

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

xu(nf) vector containing upper bounds for variables.

real(kind=wp), intent(out), dimension(nc+1) :: cf

cf(nc+1) vector containing values of the constraint functions.

integer, intent(in), dimension(nc) :: constraint_type

ic(nc) vector containing types of constraints:

  • ic(kc) = 0 -- constraint cf(kc) is not used.
  • ic(kc) = 1 -- lower constraint cl(kc) <= cf(kc).
  • ic(kc) = 2 -- upper constraint cf(kc) <= cu(kc).
  • ic(kc) = 3 -- two side constraint cl(kc) <= cf(kc) <= cu(kc).
  • ic(kc) = 5 -- equality constraint cf(kc) == cl(kc).
real(kind=wp), intent(in), dimension(nc) :: cl

cl(nc) vector containing lower bounds for constraint functions.

real(kind=wp), intent(in), dimension(nc) :: cu

cu(nc) vector containing upper bounds for constraint functions.

integer, intent(in), dimension(6) :: ipar

integer paremeters:

  • ipar(1) maximum number of iterations.
  • ipar(2) maximum number of function evaluations.
  • ipar(3) this parameter is not used in the subroutine psqp.
  • ipar(4) this parameter is not used in the subroutine psqp.
  • ipar(5) variable metric update used. ipar(5)=1 - the bfgs update. ipar(5)=2 - the hoshino update.
  • ipar(6) correction of the variable metric update if a negative curvature occurs. ipar(6)=1 - no correction. ipar(6)=2 - powell's correction.
real(kind=wp), intent(in), dimension(5) :: rpar

real parameters:

  • rpar(1) -- maximum stepsize.
  • rpar(2) -- tolerance for change of variables.
  • rpar(3) -- tolerance for constraint violations.
  • rpar(4) -- tolerance for the gradient of the lagrangian function.
  • rpar(5) -- penalty coefficient.
real(kind=wp), intent(out) :: f

value of the objective function.

real(kind=wp), intent(out) :: gmax

maximum partial derivative of the lagrangian function.

real(kind=wp), intent(out) :: cmax

maximum constraint violation.

integer, intent(in) :: iprnt

print specification:

  • iprnt=0 -- no print.
  • abs(iprnt)=1 -- print of final results.
  • abs(iprnt)=2 -- print of final results and iterations.
  • iprnt>0 -- basic final results.
  • iprnt<0 -- extended final results.
integer, intent(out) :: iterm

variable that indicates the cause of termination.

  • iterm=1 -- if abs(x-xo) was less than or equal to tolx in mtesx (usually two) subsequent iterations.
  • iterm=2 -- if abs(f-fo) was less than or equal to tolf in mtesf (usually two) subsequent iterations.
  • iterm=3 -- if f is less than or equal to tolb.
  • iterm=4 -- if gmax is less than or equal to tolg.
  • iterm=11 -- if nit exceeded mit. iterm=12 - if nfv exceeded mfv.
  • iterm=13 -- if nfg exceeded mfg. iterm<0 - if the method failed.
  • iterm=-6 -- then the termination criterion has not been satisfied, but the point obtained if usually acceptable.
procedure(obj_func) :: obj

computation of the value of the objective function

procedure(dobj_func) :: dobj

computation of the gradient of the objective function

procedure(con_func) :: con

computation of the value of the constraint function

procedure(dcon_func) :: dcon

computation of the gradient of the constraint function

procedure(report_f), optional :: report

iteration report function. Note: this is independent of iprnt. If this function is associated, each iteration will be reported


Calls

proc~~psqpn~~CallsGraph proc~psqpn psqp_class%psqpn proc~psqp psqp_class%psqp proc~psqpn->proc~psqp proc~bfgs_variable_metric_update bfgs_variable_metric_update proc~psqp->proc~bfgs_variable_metric_update proc~compute_augmented_lagrangian compute_augmented_lagrangian proc~psqp->proc~compute_augmented_lagrangian proc~compute_con_and_dcon psqp_class%compute_con_and_dcon proc~psqp->proc~compute_con_and_dcon proc~compute_new_penalty_parameters compute_new_penalty_parameters proc~psqp->proc~compute_new_penalty_parameters proc~compute_obj_and_dobj psqp_class%compute_obj_and_dobj proc~psqp->proc~compute_obj_and_dobj proc~dual_range_space_qp dual_range_space_qp proc~psqp->proc~dual_range_space_qp proc~dual_range_space_quad_prog psqp_class%dual_range_space_quad_prog proc~psqp->proc~dual_range_space_quad_prog proc~extended_line_search psqp_class%extended_line_search proc~psqp->proc~extended_line_search proc~mxdsmi mxdsmi proc~psqp->proc~mxdsmi proc~mxvcop mxvcop proc~psqp->proc~mxvcop proc~mxvdir mxvdir proc~psqp->proc~mxvdir proc~mxvdot mxvdot proc~psqp->proc~mxvdot proc~mxvina mxvina proc~psqp->proc~mxvina proc~mxvmax mxvmax proc~psqp->proc~mxvmax proc~mxvset mxvset proc~psqp->proc~mxvset proc~test_simple_bound test_simple_bound proc~psqp->proc~test_simple_bound proc~transform_incompatible_qp_subproblem transform_incompatible_qp_subproblem proc~psqp->proc~transform_incompatible_qp_subproblem proc~bfgs_variable_metric_update->proc~mxvcop proc~bfgs_variable_metric_update->proc~mxvdir proc~bfgs_variable_metric_update->proc~mxvdot proc~mxdpgb mxdpgb proc~bfgs_variable_metric_update->proc~mxdpgb proc~mxdpgp mxdpgp proc~bfgs_variable_metric_update->proc~mxdpgp proc~mxdpgs mxdpgs proc~bfgs_variable_metric_update->proc~mxdpgs proc~mxdpgu mxdpgu proc~bfgs_variable_metric_update->proc~mxdpgu proc~mxvdif mxvdif proc~bfgs_variable_metric_update->proc~mxvdif proc~mxvscl mxvscl proc~bfgs_variable_metric_update->proc~mxvscl proc~compute_con_and_dcon->proc~mxvcop proc~mxvneg mxvneg proc~compute_obj_and_dobj->proc~mxvneg proc~dual_range_space_qp->proc~mxvdir proc~dual_range_space_qp->proc~mxvdif proc~mxvsav mxvsav proc~dual_range_space_qp->proc~mxvsav proc~dual_range_space_quad_prog->proc~mxvcop proc~dual_range_space_quad_prog->proc~mxvdir proc~dual_range_space_quad_prog->proc~mxvina proc~determine_new_active_linear_constr determine_new_active_linear_constr proc~dual_range_space_quad_prog->proc~determine_new_active_linear_constr proc~determine_new_active_simple_bound determine_new_active_simple_bound proc~dual_range_space_quad_prog->proc~determine_new_active_simple_bound proc~dual_range_space_quad_prog->proc~mxdpgb proc~mxdpgf mxdpgf proc~dual_range_space_quad_prog->proc~mxdpgf proc~mxdprb mxdprb proc~dual_range_space_quad_prog->proc~mxdprb proc~mxdsmm mxdsmm proc~dual_range_space_quad_prog->proc~mxdsmm proc~mxvinv mxvinv proc~dual_range_space_quad_prog->proc~mxvinv proc~dual_range_space_quad_prog->proc~mxvneg proc~ops_after_constr_deletion psqp_class%ops_after_constr_deletion proc~dual_range_space_quad_prog->proc~ops_after_constr_deletion proc~update_tri_decomp_general update_tri_decomp_general proc~dual_range_space_quad_prog->proc~update_tri_decomp_general proc~line_search_interpolation line_search_interpolation proc~extended_line_search->proc~line_search_interpolation proc~determine_new_active_linear_constr->proc~mxvdot proc~mxdpgu->proc~mxvscl proc~update_tri_decomp_orthogonal update_tri_decomp_orthogonal proc~ops_after_constr_deletion->proc~update_tri_decomp_orthogonal proc~update_tri_decomp_general->proc~mxvcop proc~update_tri_decomp_general->proc~mxvdot proc~update_tri_decomp_general->proc~mxvset proc~update_tri_decomp_general->proc~mxdpgb proc~update_tri_decomp_general->proc~mxdprb proc~update_tri_decomp_general->proc~mxdsmm proc~mxdsmv mxdsmv proc~update_tri_decomp_general->proc~mxdsmv proc~update_tri_decomp_orthogonal->proc~mxvcop proc~update_tri_decomp_orthogonal->proc~mxvset proc~mxvort mxvort proc~update_tri_decomp_orthogonal->proc~mxvort proc~mxvrot mxvrot proc~update_tri_decomp_orthogonal->proc~mxvrot

Source Code

   subroutine psqpn(me, nf, nb, nc, x, bound_type, xl, xu, cf, constraint_type, &
                    cl, cu, ipar, rpar, f, gmax, &
                    cmax, iprnt, iterm, obj, dobj, con, dcon, report)

      class(psqp_class), intent(inout) :: me
      integer, intent(in) :: nf  !! number of variables
      integer, intent(in) :: nb  !! choice of simple bounds.
                                 !!
                                 !! * `nb=0` -- simple bounds suppressed.
                                 !! * `nb>0` -- simple bounds accepted.
      integer, intent(in) :: nc  !! number of general nonlinear constraints.
      real(wp),dimension(nf),intent(inout) :: x !! x(nf) vector of variables.
      integer, dimension(nf), intent(in) :: bound_type !! ix(nf) vector containing types of bounds.
                                                       !!
                                                       !! * `ix(i) = 0` -- variable x(i) is unbounded.
                                                       !! * `ix(i) = 1` -- lower bound xl(i) <= x(i).
                                                       !! * `ix(i) = 2` -- upper bound x(i) <= xu(i).
                                                       !! * `ix(i) = 3` -- two side bound xl(i) <= x(i) <= xu(i).
                                                       !! * `ix(i) = 5` -- variable x(i) is fixed.
      real(wp),dimension(nf),intent(in) :: xl   !! xl(nf) vector containing lower bounds for variables.
      real(wp),dimension(nf),intent(in) :: xu   !! xu(nf) vector containing upper bounds for variables.
      real(wp),dimension(nc+1),intent(out) :: cf !! cf(nc+1) vector containing values of the constraint functions.
      integer,dimension(nc),intent(in) :: constraint_type  !! ic(nc) vector containing types of constraints:
                                                           !!
                                                           !! * `ic(kc) = 0` -- constraint cf(kc) is not used.
                                                           !! * `ic(kc) = 1` -- lower constraint cl(kc) <= cf(kc).
                                                           !! * `ic(kc) = 2` -- upper constraint cf(kc) <= cu(kc).
                                                           !! * `ic(kc) = 3` -- two side constraint cl(kc) <= cf(kc) <= cu(kc).
                                                           !! * `ic(kc) = 5` -- equality constraint cf(kc) == cl(kc).
      real(wp),dimension(nc),intent(in) :: cl  !! cl(nc) vector containing lower bounds for constraint functions.
      real(wp),dimension(nc),intent(in) :: cu  !! cu(nc) vector containing upper bounds for constraint functions.
      integer,dimension(6),intent(in) :: ipar  !! integer paremeters:
                                               !!
                                               !! * `ipar(1)`  maximum number of iterations.
                                               !! * `ipar(2)`  maximum number of function evaluations.
                                               !! * `ipar(3)`  this parameter is not used in the subroutine psqp.
                                               !! * `ipar(4)`  this parameter is not used in the subroutine psqp.
                                               !! * `ipar(5)`  variable metric update used.
                                               !!   `ipar(5)=1` - the bfgs update.
                                               !!   `ipar(5)=2` - the hoshino update.
                                               !! * `ipar(6)`  correction of the variable metric update if a negative
                                               !!   curvature occurs.
                                               !!   `ipar(6)=1` - no correction.
                                               !!   `ipar(6)=2` - powell's correction.
      real(wp),dimension(5),intent(in) :: rpar   !! real parameters:
                                                 !!
                                                 !! * `rpar(1)` -- maximum stepsize.
                                                 !! * `rpar(2)` -- tolerance for change of variables.
                                                 !! * `rpar(3)` -- tolerance for constraint violations.
                                                 !! * `rpar(4)` -- tolerance for the gradient of the lagrangian function.
                                                 !! * `rpar(5)` -- penalty coefficient.
      real(wp),intent(out) :: f     !! value of the objective function.
      real(wp),intent(out) :: gmax  !! maximum partial derivative of the lagrangian function.
      real(wp),intent(out) :: cmax  !! maximum constraint violation.
      integer,intent(in) :: iprnt  !! print specification:
                                   !!
                                   !! * `iprnt=0`      -- no print.
                                   !! * `abs(iprnt)=1` -- print of final results.
                                   !! * `abs(iprnt)=2` -- print of final results and iterations.
                                   !! * `iprnt>0`      -- basic final results.
                                   !! * `iprnt<0`      -- extended final results.
      integer, intent(out) :: iterm !! variable that indicates the cause of termination.
                                    !!
                                    !! * `iterm=1`  -- if abs(x-xo) was less than or equal to tolx in mtesx (usually two) subsequent iterations.
                                    !! * `iterm=2`  -- if abs(f-fo) was less than or equal to tolf in mtesf (usually two) subsequent iterations.
                                    !! * `iterm=3`  -- if f is less than or equal to tolb.
                                    !! * `iterm=4`  -- if gmax is less than or equal to tolg.
                                    !! * `iterm=11` -- if nit exceeded mit. iterm=12 - if nfv exceeded mfv.
                                    !! * `iterm=13` -- if nfg exceeded mfg. iterm<0 - if the method failed.
                                    !! * `iterm=-6` -- then the termination criterion has not been satisfied, but the point obtained if usually acceptable.
      procedure(obj_func)  :: obj  !! computation of the value of the objective function
      procedure(dobj_func) :: dobj !! computation of the gradient of the objective function
      procedure(con_func)  :: con  !! computation of the value of the constraint function
      procedure(dcon_func) :: dcon !! computation of the gradient of the constraint function
      procedure(report_f),optional :: report !! iteration report function. Note: this is independent of `iprnt`. 
                                             !! If this function is associated, each iteration will be reported

      integer :: lcfd, lcfo, lcg, lcp, lcr, lcz, lg, lgc, lgf, lgo, lh, lia, ls, lxo
      integer, dimension(:), allocatable :: ia
      real(wp), dimension(:), allocatable :: ra
      integer,dimension(:),allocatable :: ic !! local copy of `constraint_type` since it is modified
      integer,dimension(:),allocatable :: ix !! local copy of `bound_type` since it is modified

      ! set the functions:
      me%obj  => obj
      me%dobj => dobj
      me%con  => con
      me%dcon => dcon

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

      allocate (ia(nf), ra((nf + nc + 8)*nf + 3*nc + 1))
      allocate (ic(nc))
      allocate (ix(nf))

      ic = constraint_type ! make a copy of input
      ix = bound_type ! make a copy of input
      lcg = 1
      lcfo = lcg + nf*nc
      lcfd = lcfo + nc + 1
      lgc = lcfd + nc
      lcr = lgc + nf
      lcz = lcr + nf*(nf + 1)/2
      lcp = lcz + nf
      lgf = lcp + nc
      lg = lgf + nf
      lh = lg + nf
      ls = lh + nf*(nf + 1)/2
      lxo = ls + nf
      lgo = lxo + nf
      lia = 1
      call me%psqp(nf, nb, nc, x, ix, xl, xu, cf, ic, cl, cu, ra, ra(lcfo), ra(lcfd), &
                   ra(lgc), ia, ra(lcr), ra(lcz), ra(lcp), ra(lgf), ra(lg), ra(lh), &
                   ra(ls), ra(lxo), ra(lgo), rpar(1), rpar(2), rpar(3), rpar(4), &
                   rpar(5), cmax, gmax, f, ipar(1), ipar(2), ipar(5), ipar(6), &
                   iprnt, iterm)

      deallocate (ia, ra, ic, ix)

   end subroutine psqpn