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

NAG AD Library Introduction
Example description
!   E04GB_A1W_F Example Program Text
!   Mark 28.5 Release. NAG Copyright 2022.
    Module e04gb_a1w_fe_mod

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

!     .. Use Statements ..
      Use iso_c_binding, Only: c_ptr
      Use nagad_library, Only: nagad_a1w_w_rtype, Operator (+),                &
                               Assignment (=), Operator (-), Operator (/),     &
                               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_a1w_w_rtype), Public, Save :: t(m,nt), y(m)
    Contains
      Subroutine lsqgrd(ad_handle,m,n,fvec,fjac,ldfjac,g_a1w)
!       Routine to evaluate gradient of the sum of squares

!       .. Use Statements ..
        Use nagad_library, Only: dgemv_a1w
!       .. Scalar Arguments ..
        Type (c_ptr), Intent (Inout)   :: ad_handle
        Integer, Intent (In)           :: ldfjac, m, n
!       .. Array Arguments ..
        Type (nagad_a1w_w_rtype), Intent (In) :: fjac(ldfjac,n), fvec(m)
        Type (nagad_a1w_w_rtype), Intent (Out) :: g_a1w(n)
!       .. Local Scalars ..
        Type (nagad_a1w_w_rtype)       :: alpha, beta
        Character (1), Save            :: trans = 'T'
!       .. Executable Statements ..
        alpha = 1.0_nag_wp
        beta = 0.0_nag_wp
        Call dgemv_a1w(trans,m,n,alpha,fjac,ldfjac,fvec,inc1,beta,g_a1w,inc1)

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

        Return

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

!       Routine to evaluate the residuals
!       Always use algorithmic differentiation here for simplicity

!       .. Scalar Arguments ..
        Type (c_ptr), Intent (Inout)   :: ad_handle
        Integer, Intent (Inout)        :: iflag
        Integer, Intent (In)           :: ldfjac, m, n
!       .. Array Arguments ..
        Type (nagad_a1w_w_rtype), Intent (Inout) :: fjac(ldfjac,n), ruser(*)
        Type (nagad_a1w_w_rtype), Intent (Out) :: fvec(m)
        Type (nagad_a1w_w_rtype), Intent (In) :: xc(n)
        Integer, Intent (Inout)        :: iuser(*)
!       .. Local Scalars ..
        Type (nagad_a1w_w_rtype)       :: denom, denom2
        Integer                        :: i
!       .. Executable Statements ..
!       Get computational mode
        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)
          If (iflag/=0) Then
            fjac(i,1) = 1.0_nag_wp
            denom2 = -1.0_nag_wp/(denom*denom)
            fjac(i,2) = t(i,1)*t(i,2)*denom2
            fjac(i,3) = t(i,1)*t(i,3)*denom2
          End If
        End Do
        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_a1w
!       .. 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_a1w_w_rtype), Intent (In) :: fjac(ldfjac,n), fvec(m),      &
                                          s(n), xc(n)
        Type (nagad_a1w_w_rtype), Intent (Inout) :: ruser(*)
        Integer, Intent (Inout)        :: iuser(*)
!       .. Local Scalars ..
        Type (nagad_a1w_w_rtype)       :: fsumsq, gtg
        Integer                        :: j
!       .. Local Arrays ..
        Type (nagad_a1w_w_rtype)       :: g_a1w(ndec)
!       .. Executable Statements ..
        fsumsq = ddot_a1w(m,fvec,inc1,fvec,inc1)

        Call lsqgrd(ad_handle,m,n,fvec,fjac,ldfjac,g_a1w)

        gtg = ddot_a1w(n,g_a1w,inc1,g_a1w,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_a1w(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 e04gb_a1w_fe_mod
    Program e04gb_a1w_fe

!     E04GB_A1W_F Example Main Program

!     .. Use Statements ..
      Use e04gb_a1w_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: e04gb_a1w_f, nagad_a1w_get_derivative,          &
                               nagad_a1w_inc_derivative,                       &
                               nagad_a1w_ir_create => x10za_a1w_f,             &
                               nagad_a1w_ir_interpret_adjoint_sparse,          &
                               nagad_a1w_ir_register_variable,                 &
                               nagad_a1w_ir_remove, nagad_a1w_ir_zero_adjoints &
                               , nagad_a1w_w_rtype, nagad_algorithmic,         &
                               x10aa_a1w_f, x10ab_a1w_f, x10ac_a1w_f,          &
                               x10ae_a1w_f, x10af_a1w_f, Assignment (=)
      Use nag_library, Only: nag_wp, x02ajf
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Type (c_ptr)                     :: ad_handle
      Type (nagad_a1w_w_rtype)         :: eta, fsumsq, stepmx, xtol
      Real (Kind=nag_wp)               :: yr
      Integer                          :: i, ifail, iprint, maxcal, mode, nf,  &
                                          niter, selct
!     .. Local Arrays ..
      Type (nagad_a1w_w_rtype)         :: fjac(m,n), fvec(m), g_a1w(n),        &
                                          ruser(1), s(n), v(ldv,n), x(n)
      Real (Kind=nag_wp)               :: dfdy(m), dxdy(m), tr(nt)
      Integer                          :: iuser(1)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: real, sqrt
!     .. Executable Statements ..
      Write (nout,*) 'E04GB_A1W_F Example Program Results'

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

!     Read the computational mode: 1 = algorithmic, 2 = symbolic
      Read (nin,*) mode

!     Observations of T_j (j = 1: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 = 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
      Do i = 1, n
        x(i) = 0.5_nag_wp*real(i,kind=nag_wp)
      End Do

!     Create AD tape
      Call nagad_a1w_ir_create

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

!     Set computational mode
      ifail = 0
      Call x10ac_a1w_f(ad_handle,mode,ifail)
      Call x10ae_a1w_f(ad_handle,-1,ifail)
      Call x10af_a1w_f(ad_handle,-1,ifail)

!     Register variables to differentiate w.r.t.
      Call nagad_a1w_ir_register_variable(y)
      Call nagad_a1w_ir_register_variable(t)

!     Solve the problem
      selct = 2
      ifail = -1

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

      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(ad_handle,m,n,fvec,fjac,ldfjac,g_a1w)

        Write (nout,99998) 'Estim gradien  = ', g_a1w(1:n)%value
        Write (nout,*) '                    (machine dependent)'
        Write (nout,*) 'Residuals:'
        Write (nout,99997) fvec(1:m)%value
      Case (:-1)
        Write (nout,*)
        Write (nout,99995) 'Routine e04gb_a1w_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))

      Write (nout,*)
      Write (nout,*) ' Derivatives calculated: First order adjoints'
      If (mode==nagad_algorithmic) Then
        Write (nout,*) ' Computational mode    : algorithmic'
      Else
        Write (nout,*) ' Computational mode    : symbolic'
      End If
      Write (nout,*)
      Write (nout,*) ' Derivatives:'

!     Setup evaluation of derivatives via adjoints
      Call nagad_a1w_inc_derivative(fsumsq,1.0_nag_wp)
      ifail = 0
      Call nagad_a1w_ir_interpret_adjoint_sparse(ifail)

!     Get derivatives
      dfdy = nagad_a1w_get_derivative(y)

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

      Do i = 1, n
        Call nagad_a1w_ir_zero_adjoints
        Call nagad_a1w_inc_derivative(x(i),1.0_nag_wp)
        ifail = 0
        Call nagad_a1w_ir_interpret_adjoint_sparse(ifail)

!       Get derivatives
        dxdy = nagad_a1w_get_derivative(y)

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

100   Continue
!     Remove computational data object and tape
      Call x10ab_a1w_f(ad_handle,ifail)
      Call nagad_a1w_ir_remove

    End Program e04gb_a1w_fe