NAG Library Manual, Mark 29.3
Interfaces:  FL   CL   CPP   AD 

NAG AD Library Introduction
Example description
!   E04FC_T1W_F Example Program Text
!   Mark 29.3 Release. NAG Copyright 2023.

    Module e04fc_t1w_fe_mod

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

!     .. Use Statements ..
      Use iso_c_binding, Only: c_ptr
      Use nagad_library, Only: nagad_t1w_w_rtype, Operator (-), Operator (+),  &
                               Operator (/), Assignment (=), Operator (*)
      Use nag_library, Only: nag_wp
!     .. Implicit None Statement ..
      Implicit None
!     .. Accessibility Statements ..
      Private
      Public                           :: lsqfun, lsqgrd, lsqmon
!     .. Parameters ..
      Integer, Parameter               :: inc1 = 1
      Integer, Parameter, Public       :: m = 15, n = 3, nin = 5, nout = 6,    &
                                          nt = 3
      Integer, Parameter, Public       :: ldfjac = m
      Integer, Parameter, Public       :: ldv = n
!     .. Local Arrays ..
      Type (nagad_t1w_w_rtype), 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 nagad_library, Only: dgemv_t1w
!       .. Scalar Arguments ..
        Integer, Intent (In)           :: ldfjac, m, n
!       .. Array Arguments ..
        Type (nagad_t1w_w_rtype), Intent (In) :: fjac(ldfjac,n), fvec(m)
        Type (nagad_t1w_w_rtype), Intent (Out) :: g(n)
!       .. Local Scalars ..
        Type (nagad_t1w_w_rtype)       :: alpha, beta
        Character (1), Save            :: trans = 'T'
!       .. Executable Statements ..
        alpha = 1.0_nag_wp
        beta = 0.0_nag_wp
        Call dgemv_t1w(trans,m,n,alpha,fjac,ldfjac,fvec,inc1,beta,g,inc1)

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

        Return

      End Subroutine lsqgrd
      Subroutine lsqfun(ad_handle,iflag,m,n,xc,fvec,iuser,ruser)

!       Routine to evaluate the residuals

!       .. Scalar Arguments ..
        Type (c_ptr), Intent (Inout)   :: ad_handle
        Integer, Intent (Inout)        :: iflag
        Integer, Intent (In)           :: m, n
!       .. Array Arguments ..
        Type (nagad_t1w_w_rtype), Intent (Out) :: fvec(m)
        Type (nagad_t1w_w_rtype), Intent (Inout) :: ruser(*)
        Type (nagad_t1w_w_rtype), Intent (In) :: xc(n)
        Integer, Intent (Inout)        :: iuser(*)
!       .. Executable Statements ..
        fvec(1:m) = xc(1) + t(1:m,1)/(xc(2)*t(1:m,2)+xc(3)*t(1:m,3)) - y(1:m)

        Return

      End Subroutine lsqfun
      Subroutine lsqmon(ad_handle,m,n,xc,fvec,fjac,ldfjac,s,igrade,niter,nf,   &
        iuser,ruser)
!       Monitoring routine

!       .. Use Statements ..
        Use nagad_library, Only: ddot_t1w
!       .. Parameters ..
        Integer, Parameter             :: ndec = 3
!       .. Scalar Arguments ..
        Type (c_ptr), Intent (Inout)   :: ad_handle
        Integer, Intent (In)           :: igrade, ldfjac, m, n, nf, niter
!       .. Array Arguments ..
        Type (nagad_t1w_w_rtype), Intent (In) :: fjac(ldfjac,n), fvec(m),      &
                                          s(n), xc(n)
        Type (nagad_t1w_w_rtype), Intent (Inout) :: ruser(*)
        Integer, Intent (Inout)        :: iuser(*)
!       .. Local Scalars ..
        Type (nagad_t1w_w_rtype)       :: fsumsq, gtg
        Integer                        :: j
!       .. Local Arrays ..
        Type (nagad_t1w_w_rtype)       :: g(ndec)
!       .. Executable Statements ..
!       The NAG name equivalent of ddot_t1w is f06ea_t1w_f
        fsumsq = ddot_t1w(m,fvec,inc1,fvec,inc1)

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

        gtg = ddot_t1w(n,g,inc1,g,inc1)

        Write (nout,*)
        Write (nout,*)                                                         &
          '  Itn      F evals        SUMSQ             GTG        Grade'
        Write (nout,99999) niter, nf, fsumsq%value, gtg%value, igrade
        Write (nout,*)
        Write (nout,*)                                                         &
          '       x                    g           Singular values'

        Write (nout,99998)(xc(j)%value,g(j)%value,s(j)%value,j=1,n)

        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 e04fc_t1w_fe_mod
    Program e04fc_t1w_fe

!     E04FC_T1W_F Example Main Program

!     .. Use Statements ..
      Use e04fc_t1w_fe_mod, Only: ldfjac, ldv, lsqfun, lsqgrd, lsqmon, m, n,   &
                                  nin, nout, nt, t, y
      Use iso_c_binding, Only: c_ptr
      Use nagad_library, Only: e04fc_t1w_f, nagad_t1w_set_derivative,          &
                               nagad_t1w_w_rtype, x10aa_t1w_f, x10ab_t1w_f,    &
                               Assignment (=)
      Use nag_library, Only: nag_wp, x02ajf
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Type (c_ptr)                     :: ad_handle
      Type (nagad_t1w_w_rtype)         :: eta, fsumsq, stepmx, xtol
      Real (Kind=nag_wp)               :: yr
      Integer                          :: i, ifail, iprint, j, maxcal, nf,     &
                                          niter
!     .. Local Arrays ..
      Type (nagad_t1w_w_rtype)         :: fjac(m,n), fvec(m), g(n), ruser(1),  &
                                          s(n), v(ldv,n), x(n)
      Real (Kind=nag_wp)               :: dfdy(m), dxdy(n,m), tr(nt)
      Integer                          :: iuser(1)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: real, sqrt
!     .. Executable Statements ..
      Write (nout,*) 'E04FC_T1W_F Example Program Results'

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

!     Observations of t_j (j = 1, 2, ..., nt) are held in t(1:m, J)
      Do i = 1, m
        Read (nin,*) yr, tr(1:nt)
        y(i) = yr
        t(i,1:nt) = tr(1:nt)
      End Do

!     Set IPRINT to 1 to obtain output from LSQMON at each iteration
      iprint = -1
      maxcal = 400*n
      eta = 0.5_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

!     Create AD configuration data object
      ifail = 0
      Call x10aa_t1w_f(ad_handle,ifail)

      Do i = 1, m
        Call nagad_t1w_set_derivative(y(i),1.0_nag_wp)
!       Solve the problem
!       Set up the starting point
        Do j = 1, n
          x(j) = 0.5_nag_wp*real(j,kind=nag_wp)
        End Do

        ifail = -1
        Call e04fc_t1w_f(ad_handle,m,n,lsqfun,lsqmon,iprint,maxcal,eta,xtol,   &
          stepmx,x,fsumsq,fvec,fjac,ldfjac,s,v,ldv,niter,nf,iuser,ruser,ifail)

        y(i)%tangent = 0.0_nag_wp
!       Get derivatives
        dfdy(i) = fsumsq%tangent
        dxdy(1:n,i) = x(1:n)%tangent
      End Do

!     Primal results
      Select Case (ifail)
      Case (0,2:)
        Write (nout,*)
        Write (nout,99999) 'Sum of squares = ', fsumsq%value
        Write (nout,99999) 'Solution point = ', x(1:n)%value

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

        Write (nout,99998) 'Estim gradient = ', g(1:n)%value
        Write (nout,*) '                    (machine dependent)'
        Write (nout,*) 'Residuals:'
        Write (nout,99997) fvec(1:m)%value
      Case (:-1)
        Write (nout,99996) 'Routine e04fc_t1w_f failed with ifail = ', ifail
        Go To 100
      End Select

99999 Format (1X,A,3F12.4)
99998 Format (1X,A,1P,3E12.3)
99997 Format (3(1X,1P,E11.1))
99996 Format (/,1X,A,I0)

      Write (nout,*)
      Write (nout,*) ' Derivatives calculated: First order tangents'
      Write (nout,*) ' Computational mode    : algorithmic'

      Write (nout,*)
      Write (nout,*) ' Derivatives:'
      Write (nout,*)
      Write (nout,*) '  sum of squares w.r.t y:'
      Write (nout,99995) dfdy
99995 Format (3(1X,1P,E11.1))

      Do i = 1, n
        Write (nout,*)
        Write (nout,99994) '  dx(', i, ')/dy:'
        Write (nout,99995) dxdy(i,1:m)
      End Do
99994 Format (1X,A,I0,A)

100   Continue
!     Remove computational data object
      Call x10ab_t1w_f(ad_handle,ifail)

    End Program e04fc_t1w_fe