!   Mark 25 Release. NAG Copyright 2014.

    Module e04mxfe_mod

!     E04MXF Example Program Module:
!            Parameters and User-defined Routines

!     .. Use Statements ..
      Use nag_library, Only: nag_wp
!     .. Implicit None Statement ..
      Implicit None
!     .. Accessibility Statements ..
      Private
      Public                               :: qphx
    Contains
      Subroutine qphx(ncolh,x,hx,nstate,cuser,iuser,ruser)

!       Subroutine to compute H*x. 

!       Note: IUSER and RUSER contain the following data:
!       RUSER(1:NNZH) = H(1:NNZH)
!       IUSER(1:NCOLH+1) = ICCOLH(1:NCOLH+1)
!       IUSER(NCOLH+2:NNZH+NCOLH+1) = IROWH(1:NNZH)

!       .. Scalar Arguments ..
        Integer, Intent (In)                 :: ncolh, nstate
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out)     :: hx(ncolh)
        Real (Kind=nag_wp), Intent (Inout)   :: ruser(*)
        Real (Kind=nag_wp), Intent (In)      :: x(ncolh)
        Integer, Intent (Inout)              :: iuser(*)
        Character (8), Intent (Inout)        :: cuser(*)
!       .. Local Scalars ..
        Integer                              :: end, icol, idx, irow, start
!       .. Executable Statements ..
        hx(1:ncolh) = 0.0E0_nag_wp

        Do icol = 1, ncolh
          start = iuser(icol)
          end = iuser(icol+1) - 1

          Do idx = start, end
            irow = iuser(ncolh+1+idx)
            hx(irow) = hx(irow) + x(icol)*ruser(idx)
            If (irow/=icol) Then
              hx(icol) = hx(icol) + x(irow)*ruser(idx)
            End If

          End Do

        End Do

        Return
      End Subroutine qphx
    End Module e04mxfe_mod

    Program e04mxfe

!     .. Use Statements ..
      Use nag_library, Only: e04mxf, e04npf, e04nqf, e04nsf, e04ntf, nag_wp,   &
                             x04acf, x04adf
      Use e04mxfe_mod, Only: qphx
!     .. Implicit None Statement ..
      Implicit None
!     .. Parameters ..
      Integer, Parameter                   :: lencw = 600, leniw = 600,        &
                                              lenrw = 600, mpslst = 1,         &
                                              nin = 7, nout = 6
      Logical, Parameter                   :: readints = .False.
      Character (*), Parameter             :: fname = 'e04mxfe.opt'
!     .. Local Scalars ..
      Real (Kind=nag_wp)                   :: obj, objadd, sinf
      Integer                              :: i, ifail, iobj, lenc, lintvar,   &
                                              m, maxlintvar, maxm, maxn,       &
                                              maxncolh, maxnnz, maxnnzh,       &
                                              minmax, mode, n, ncolh, ninf,    &
                                              nname, nnz, nnzh, ns
      Character (1)                        :: start
!     .. Local Arrays ..
      Real (Kind=nag_wp), Allocatable      :: a(:), bl(:), bu(:), c(:), h(:),  &
                                              pi(:), rc(:), ruser(:), rw(:),   &
                                              x(:)
      Integer, Allocatable                 :: helast(:), hs(:), iccola(:),     &
                                              iccolh(:), intvar(:), irowa(:),  &
                                              irowh(:), iuser(:), iw(:)
      Character (8), Allocatable           :: crname(:), cw(:)
      Character (8)                        :: cuser(1), pnames(5)
!     .. Intrinsic Procedures ..
      Intrinsic                            :: max, min
!     .. Executable Statements ..
      Write (nout,*) 'E04MXF Example Program Results'
      Flush (nout)

!     Initialize
      pnames(1:5) = '        '
      maxm = 0
      maxn = 0
      maxnnz = 0
      maxnnzh = 0
      maxncolh = 0
      maxlintvar = 0

!     Open the data file for reading
      mode = 0
      ifail = 0
      Call x04acf(nin,fname,mode,ifail)

!     Call e04mxf in query mode
      Allocate (a(maxnnz),irowa(maxnnz),iccola(maxn+1),bl(maxn+maxm), &
        bu(maxn+maxm),crname(maxn+maxm),h(maxnnzh),irowh(maxnnzh), &
        iccolh(maxncolh+1),intvar(maxlintvar))
      ifail = 0
      Call e04mxf(nin,maxn,maxm,maxnnz,maxncolh,maxnnzh,maxlintvar,mpslst,n,m, &
        nnz,ncolh,nnzh,lintvar,iobj,a,irowa,iccola,bl,bu,pnames,nname,crname, &
        h,irowh,iccolh,minmax,intvar,ifail)
      Deallocate (a,irowa,iccola,bl,bu,crname,h,irowh,iccolh,intvar)

!     Close the data file
      ifail = 0
      Call x04adf(nin,ifail)

!     set maxm maxn and maxnnz 
      maxm = m
      maxn = n
      maxnnz = nnz
      maxnnzh = nnzh
      maxncolh = ncolh
      If (readints) Then
        maxlintvar = lintvar
      Else
        maxlintvar = -1
      End If

!     Allocate memory
      Allocate (irowa(maxnnz),iccola(maxn+1),a(maxnnz),bl(maxn+maxm), &
        bu(maxn+maxm),crname(maxn+maxm),irowh(maxnnzh),iccolh(maxncolh+1), &
        h(maxnnzh),intvar(maxlintvar))

!     Open the data file for reading
      mode = 0
      ifail = 0
      Call x04acf(nin,fname,mode,ifail)

!     Call e04mxf to read the problem
      ifail = 0
      Call e04mxf(nin,maxn,maxm,maxnnz,maxncolh,maxnnzh,maxlintvar,mpslst,n,m, &
        nnz,ncolh,nnzh,lintvar,iobj,a,irowa,iccola,bl,bu,pnames,nname,crname, &
        h,irowh,iccolh,minmax,intvar,ifail)

!     Close the data file
      ifail = 0
      Call x04adf(nin,ifail)

!     Data has been read. Set up and run the solver

      Allocate (iw(leniw),rw(lenrw),cw(lencw))

!     Call e04npf to initialize workspace
      ifail = 0
      Call e04npf(cw,lencw,iw,leniw,rw,lenrw,ifail)

!     Call option setter e04nsf to change the direction of optimization.
!     Minimization is assumed by default.
      If (minmax==1) Then
        ifail = 0
        Call e04nsf('Maximize',cw,iw,rw,ifail)
      Else If (minmax==0) Then
        ifail = 0
        Call e04nsf('Feasible Point',cw,iw,rw,ifail)
      End If

!     By default E04NQF does not print monitoring
!     information. Set the print file unit or the summary
!     file unit to get information.
      ifail = 0
      Call e04ntf('Print file',nout,cw,iw,rw,ifail)


!     We have no explicit objective vector so set LENC = 0; the
!     objective vector is stored in row IOBJ of ACOL.
      lenc = 0
      objadd = 0.0E0_nag_wp
      start = 'C'

      Allocate (c(max(1,lenc)),helast(n+m),x(n+m),pi(m),rc(n+m),hs(n+m),iuser( &
        ncolh+1+nnzh),ruser(nnzh))

      helast(1:n+m) = 0
      hs(1:n+m) = 0
      Do i = 1, n + m
        x(i) = min(max(0.0E0_nag_wp,bl(i)),bu(i))
      End Do

      If (ncolh>0) Then
!       Store the non zeros of H in ruser for use by qphx
        ruser(1:nnzh) = h(1:nnzh)

!       Store iccolh and irowh in iuser for use by qphx
        iuser(1:ncolh+1) = iccolh(1:ncolh+1)
        iuser(ncolh+2:nnzh+ncolh+1) = irowh(1:nnzh)
      End If

!     Call e04nqf to solve the problem
      ifail = 0
      Call e04nqf(start,qphx,m,n,nnz,nname,lenc,ncolh,iobj,objadd,pnames(1),a, &
        irowa,iccola,bl,bu,c,crname,helast,hs,x,pi,rc,ns,ninf,sinf,obj,cw, &
        lencw,iw,leniw,rw,lenrw,cuser,iuser,ruser,ifail)

    End Program e04mxfe