Example description
!   E04HEF Example Program Text
!   Mark 26.2 Release. NAG Copyright 2017.
    Module e04hefe_mod

!     E04HEF 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                           :: lsqfun, lsqgrd, lsqhes, lsqmon
!     .. Parameters ..
      Integer, Parameter, Public       :: liw = 1, m = 15, n = 3, nin = 5,     &
                                          nout = 6, nt = 3
      Integer, Parameter, Public       :: lb = n*(n+1)/2
      Integer, Parameter, Public       :: ldfjac = m
      Integer, Parameter, Public       :: ldv = n
      Integer, Parameter, Public       :: lw = 7*n + m*n + 2*m + n*n
!     .. Local Arrays ..
      Real (Kind=nag_wp), Public, Save :: t(m,nt), y(m)
    Contains
      Subroutine lsqgrd(m,n,fvec,fjac,ldfjac,g)
!       Routine to evaluate gradient of the sum of squares

!       .. Use Statements ..
        Use nag_library, Only: dgemv
!       .. Scalar Arguments ..
        Integer, Intent (In)           :: ldfjac, m, n
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (In) :: fjac(ldfjac,n), fvec(m)
        Real (Kind=nag_wp), Intent (Out) :: g(n)
!       .. Executable Statements ..
!       The NAG name equivalent of dgemv is f06paf
        Call dgemv('T',m,n,1.0_nag_wp,fjac,ldfjac,fvec,1,0.0_nag_wp,g,1)

        g(1:n) = 2.0_nag_wp*g(1:n)

        Return

      End Subroutine lsqgrd
      Subroutine lsqfun(iflag,m,n,xc,fvec,fjac,ldfjac,iw,liw,w,lw)

!       Routine to evaluate the residuals and their 1st derivatives

!       .. Scalar Arguments ..
        Integer, Intent (Inout)        :: iflag
        Integer, Intent (In)           :: ldfjac, liw, lw, m, n
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Inout) :: fjac(ldfjac,n), w(lw)
        Real (Kind=nag_wp), Intent (Out) :: fvec(m)
        Real (Kind=nag_wp), Intent (In) :: xc(n)
        Integer, Intent (Inout)        :: iw(liw)
!       .. Local Scalars ..
        Real (Kind=nag_wp)             :: denom, dummy
        Integer                        :: i
!       .. Executable Statements ..
        Do i = 1, m
          denom = xc(2)*t(i,2) + xc(3)*t(i,3)
          fvec(i) = xc(1) + t(i,1)/denom - y(i)
          fjac(i,1) = 1.0_nag_wp
          dummy = -1.0_nag_wp/(denom*denom)
          fjac(i,2) = t(i,1)*t(i,2)*dummy
          fjac(i,3) = t(i,1)*t(i,3)*dummy
        End Do

        Return

      End Subroutine lsqfun
      Subroutine lsqhes(iflag,m,n,fvec,xc,b,lb,iw,liw,w,lw)

!       Routine to compute the lower triangle of the matrix B
!       (stored by rows in the array B)

!       .. Scalar Arguments ..
        Integer, Intent (Inout)        :: iflag
        Integer, Intent (In)           :: lb, liw, lw, m, n
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out) :: b(lb)
        Real (Kind=nag_wp), Intent (In) :: fvec(m), xc(n)
        Real (Kind=nag_wp), Intent (Inout) :: w(lw)
        Integer, Intent (Inout)        :: iw(liw)
!       .. Local Scalars ..
        Real (Kind=nag_wp)             :: dummy, sum22, sum32, sum33
        Integer                        :: i
!       .. Executable Statements ..
        b(1) = 0.0_nag_wp
        b(2) = 0.0_nag_wp
        sum22 = 0.0_nag_wp
        sum32 = 0.0_nag_wp
        sum33 = 0.0_nag_wp

        Do i = 1, m
          dummy = 2.0_nag_wp*t(i,1)/(xc(2)*t(i,2)+xc(3)*t(i,3))**3
          sum22 = sum22 + fvec(i)*dummy*t(i,2)**2
          sum32 = sum32 + fvec(i)*dummy*t(i,2)*t(i,3)
          sum33 = sum33 + fvec(i)*dummy*t(i,3)**2
        End Do

        b(3) = sum22
        b(4) = 0.0_nag_wp
        b(5) = sum32
        b(6) = sum33

        Return

      End Subroutine lsqhes
      Subroutine lsqmon(m,n,xc,fvec,fjac,ldfjac,s,igrade,niter,nf,iw,liw,w,lw)
!       Monitoring routine

!       .. Use Statements ..
        Use nag_library, Only: f06eaf
!       .. Parameters ..
        Integer, Parameter             :: ndec = 3
!       .. Scalar Arguments ..
        Integer, Intent (In)           :: igrade, ldfjac, liw, lw, m, n, nf,   &
                                          niter
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (In) :: fjac(ldfjac,n), fvec(m), s(n),      &
                                          xc(n)
        Real (Kind=nag_wp), Intent (Inout) :: w(lw)
        Integer, Intent (Inout)        :: iw(liw)
!       .. Local Scalars ..
        Real (Kind=nag_wp)             :: fsumsq, gtg
        Integer                        :: j
!       .. Local Arrays ..
        Real (Kind=nag_wp)             :: g(ndec)
!       .. Executable Statements ..
        fsumsq = f06eaf(m,fvec,1,fvec,1)

        Call lsqgrd(m,n,fvec,fjac,ldfjac,g)

        gtg = f06eaf(n,g,1,g,1)

        Write (nout,*)
        Write (nout,*)                                                         &
          ' Itns    F evals          SUMSQ             GTG        grade'
        Write (nout,99999) niter, nf, fsumsq, gtg, igrade
        Write (nout,*)
        Write (nout,*)                                                         &
          '       X                    G           Singular values'

        Do j = 1, n
          Write (nout,99998) xc(j), g(j), s(j)
        End Do

        Return

99999   Format (1X,I4,6X,I5,6X,1P,E13.5,6X,1P,E9.1,6X,I3)
99998   Format (1X,1P,E13.5,10X,1P,E9.1,10X,1P,E9.1)
      End Subroutine lsqmon
    End Module e04hefe_mod
    Program e04hefe

!     E04HEF Example Main Program

!     .. Use Statements ..
      Use e04hefe_mod, Only: lb, ldfjac, ldv, liw, lsqfun, lsqgrd, lsqhes,     &
                             lsqmon, lw, m, n, nin, nout, nt, t, y
      Use nag_library, Only: e04hef, e04yaf, e04ybf, nag_wp, x02ajf
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Real (Kind=nag_wp)               :: eta, fsumsq, stepmx, xtol
      Integer                          :: i, ifail, iprint, maxcal, nf, niter
!     .. Local Arrays ..
      Real (Kind=nag_wp)               :: b(lb), fjac(ldfjac,n), fvec(m),      &
                                          g(n), s(n), v(ldv,n), w(lw), x(n)
      Integer                          :: iw(liw)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: sqrt
!     .. Executable Statements ..
      Write (nout,*) 'E04HEF Example Program Results'

!     Skip heading in data file
      Read (nin,*)

!     Observations of TJ (J = 1, 2, ..., nt) are held in T(I, J)
!     (I = 1, 2, ..., m)

      Do i = 1, m
        Read (nin,*) y(i), t(i,1:nt)
      End Do

!     Set up an arbitrary point at which to check the derivatives

      x(1:nt) = (/0.19_nag_wp,-1.34_nag_wp,0.88_nag_wp/)

!     Check the 1st derivatives

      ifail = 0
      Call e04yaf(m,n,lsqfun,x,fvec,fjac,ldfjac,iw,liw,w,lw,ifail)

!     Check the evaluation of B

      ifail = 0
      Call e04ybf(m,n,lsqfun,lsqhes,x,fvec,fjac,ldfjac,b,lb,iw,liw,w,lw,ifail)

!     Continue setting parameters for E04HEF

!     Set IPRINT to 1 to obtain output from LSQMON at each iteration

      iprint = -1

      maxcal = 50*n
      eta = 0.9_nag_wp
      xtol = 10.0_nag_wp*sqrt(x02ajf())

!     We estimate that the minimum will be within 10 units of the
!     starting point

      stepmx = 10.0_nag_wp

!     Set up the starting point

      x(1:nt) = (/0.5_nag_wp,1.0_nag_wp,1.5_nag_wp/)

      ifail = -1
      Call e04hef(m,n,lsqfun,lsqhes,lsqmon,iprint,maxcal,eta,xtol,stepmx,x,    &
        fsumsq,fvec,fjac,ldfjac,s,v,ldv,niter,nf,iw,liw,w,lw,ifail)

      Select Case (ifail)
      Case (0,2:)
        Write (nout,*)
        Write (nout,99999) 'On exit, the sum of squares is', fsumsq
        Write (nout,99999) 'at the point', x(1:n)

        Call lsqgrd(m,n,fvec,fjac,ldfjac,g)

        Write (nout,99998) 'The corresponding gradient is', g(1:n)
        Write (nout,*) '                           (machine dependent)'
        Write (nout,*) 'and the residuals are'
        Write (nout,99997) fvec(1:m)
      End Select

99999 Format (1X,A,3F12.4)
99998 Format (1X,A,1P,3E12.3)
99997 Format (1X,1P,E9.1)
    End Program e04hefe