dual_range_space_qp Subroutine

private subroutine dual_range_space_qp(nf, n, x, xo, ica, cg, cz, g, go, r, f, fo, p, po, cmax, cmaxo, dmax, kd, ld, iters)

dual range space quadratic programming method for minimax approximation.

Note

This routine was formerly called pytrnd.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nf

declared number of variables.

integer, intent(inout) :: n

actual number of variables.

real(kind=wp) :: x(*)

x(nf) vector of variables.

real(kind=wp) :: xo(*)

xo(nf) saved vector of variables.

integer, intent(in) :: ica(*)

ica(nf) vector containing indices of active constraints.

real(kind=wp) :: cg(*)

cg(nf*nc) matrix whose columns are normals of the linear constraints.

real(kind=wp) :: cz(*)

cz(nf) vector of lagrange multipliers.

real(kind=wp) :: g(*)

g(nf) gradient of the lagrangian function.

real(kind=wp) :: go(*)

go(nf) saved gradient of the lagrangian function.

real(kind=wp) :: r

value of the stepsize parameter.

real(kind=wp) :: f

new value of the objective function.

real(kind=wp) :: fo

old value of the objective function.

real(kind=wp) :: p

new value of the directional derivative.

real(kind=wp) :: po

old value of the directional derivative.

real(kind=wp) :: cmax

value of the constraint violation.

real(kind=wp) :: cmaxo

saved value of the constraint violation.

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

maximum relative difference of variables.

integer :: kd
integer :: ld
integer :: iters

termination indicator for steplength determination. iters=0 for zero step.


Calls

proc~~dual_range_space_qp~~CallsGraph proc~dual_range_space_qp dual_range_space_qp proc~mxvdif mxvdif proc~dual_range_space_qp->proc~mxvdif proc~mxvdir mxvdir proc~dual_range_space_qp->proc~mxvdir proc~mxvsav mxvsav proc~dual_range_space_qp->proc~mxvsav

Called by

proc~~dual_range_space_qp~~CalledByGraph proc~dual_range_space_qp dual_range_space_qp proc~psqp psqp_class%psqp proc~psqp->proc~dual_range_space_qp proc~psqpn psqp_class%psqpn proc~psqpn->proc~psqp

Source Code

   subroutine dual_range_space_qp(nf, n, x, xo, ica, cg, cz, g, go, r, f, fo, &
                                  p, po, cmax, cmaxo, dmax, kd, ld, iters)

      integer,intent(in) :: nf  !! declared number of variables.
      integer,intent(inout) :: n  !! actual number of variables.
      integer,intent(in) :: ica(*)  !! ica(nf)  vector containing indices of active constraints.
      real(wp) :: x(*)  !! x(nf)  vector of variables.
      real(wp) :: xo(*)  !! xo(nf)  saved vector of variables.
      real(wp) :: cg(*)  !! cg(nf*nc)  matrix whose columns are normals of the linear constraints.
      real(wp) :: cz(*)  !! cz(nf)  vector of lagrange multipliers.
      real(wp) :: g(*)  !! g(nf)  gradient of the lagrangian function.
      real(wp) :: go(*)  !! go(nf)  saved gradient of the lagrangian function.
      real(wp) :: r  !! value of the stepsize parameter.
      real(wp) :: f  !! new value of the objective function.
      real(wp) :: fo  !! old value of the objective function.
      real(wp) :: p  !! new value of the directional derivative.
      real(wp) :: po  !! old value of the directional derivative.
      real(wp) :: cmax  !! value of the constraint violation.
      real(wp) :: cmaxo  !! saved value of the constraint violation.
      real(wp),intent(out) :: dmax  !! maximum relative difference of variables.
      integer :: kd  !!
      integer :: ld  !!
      integer :: iters  !! termination indicator for steplength determination.
                        !! iters=0 for zero step.

      integer :: i, j, l

      do j = 1, nf - n
         l = ica(j)
         if (l > 0) then
            call mxvdir(nf, -cz(j), cg((l - 1)*nf + 1), g, g)
         else
            l = -l
            g(l) = g(l) - cz(j)
         end if
      end do
      if (iters > 0) then
         call mxvdif(nf, x, xo, xo)
         call mxvdif(nf, g, go, go)
         po = r*po
         p = r*p
      else
         f = fo
         p = po
         cmax = cmaxo
         call mxvsav(nf, x, xo)
         call mxvsav(nf, g, go)
         ld = kd
      end if
      dmax = 0.0_wp
      do i = 1, nf
         dmax = max(dmax, abs(xo(i))/max(abs(x(i)), 1.0_wp))
      end do
      n = nf
   end subroutine dual_range_space_qp