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

NAG AD Library Introduction
Example description
!   E04WB_A1W_F Example Program Text
!   Mark 28.7 Release. NAG Copyright 2022.
    Module e04wb_a1w_fe_mod

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

!     .. Use Statements ..
      Use iso_c_binding, Only: c_ptr
      Use nagad_library, Only: exp, nagad_a1w_w_rtype, Assignment (=),         &
                               Operator (+), Operator (-), Operator (*)
      Use nag_library, Only: nag_wp
!     .. Implicit None Statement ..
      Implicit None
!     .. Accessibility Statements ..
      Private
      Public                           :: confun, objfun
!     .. Parameters ..
      Integer, Parameter, Public       :: nin = 5, nout = 6
    Contains
      Subroutine objfun(ad_handle,mode,m,n,ldfj,needfi,x,f,fjac,nstate,iuser,  &
        ruser)
!       Routine to evaluate the subfunctions and their 1st derivatives.
!       .. Scalar Arguments ..
        Type (c_ptr), Intent (Inout)   :: ad_handle
        Integer, Intent (In)           :: ldfj, m, n, needfi, nstate
        Integer, Intent (Inout)        :: mode
!       .. Array Arguments ..
        Type (nagad_a1w_w_rtype), Intent (Out) :: f(m)
        Type (nagad_a1w_w_rtype), Intent (Inout) :: fjac(ldfj,n), ruser(*)
        Type (nagad_a1w_w_rtype), Intent (In) :: x(n)
        Integer, Intent (Inout)        :: iuser(*)
!       .. Local Scalars ..
        Type (nagad_a1w_w_rtype)       :: ai, temp, x1, x2
        Integer                        :: i
        Logical                        :: mode02, mode12
!       .. Executable Statements ..
        x1 = x(1)
        x2 = x(2)

        If (mode==0 .And. needfi>0) Then
          f(needfi) = x1 + (0.49_nag_wp-x1)*exp(-x2*(ruser(needfi)-8.0_nag_wp) &
            )
        Else

          mode02 = (mode==0 .Or. mode==2)
          mode12 = (mode==1 .Or. mode==2)

          Do i = 1, m

            ai = ruser(i) - 8.0_nag_wp
            temp = exp(-x2*ai)

            If (mode02) Then
              f(i) = x1 + (0.49_nag_wp-x1)*temp
            End If

            If (mode12) Then
              fjac(i,1) = 1.0_nag_wp - temp
              fjac(i,2) = -(0.49_nag_wp-x1)*ai*temp
            End If

          End Do
        End If

        Return

      End Subroutine objfun
      Subroutine confun(ad_handle,mode,ncnln,n,ldcj,needc,x,c,cjac,nstate,     &
        iuser,ruser)
!       Routine to evaluate the nonlinear constraint and its 1st
!       derivatives.

!       .. Scalar Arguments ..
        Type (c_ptr), Intent (Inout)   :: ad_handle
        Integer, Intent (In)           :: ldcj, n, ncnln, nstate
        Integer, Intent (Inout)        :: mode
!       .. Array Arguments ..
        Type (nagad_a1w_w_rtype), Intent (Out) :: c(ncnln)
        Type (nagad_a1w_w_rtype), Intent (Inout) :: cjac(ldcj,n), ruser(*)
        Type (nagad_a1w_w_rtype), Intent (In) :: x(n)
        Integer, Intent (Inout)        :: iuser(*)
        Integer, Intent (In)           :: needc(ncnln)
!       .. Executable Statements ..
        If (nstate==1) Then

!         First call to CONFUN.  Set all Jacobian elements to zero.
!         Note that this will only work when 'Derivative Level = 3'
!         (the default; see Section 11.2).

          cjac(1:ncnln,1:n) = 0.0_nag_wp
        End If

        If (needc(1)>0) Then

          If (mode==0 .Or. mode==2) Then
            c(1) = -0.09_nag_wp - x(1)*x(2) + 0.49_nag_wp*x(2)
          End If

          If (mode==1 .Or. mode==2) Then
            cjac(1,1) = -x(2)

            cjac(1,2) = -x(1) + 0.49_nag_wp
          End If

        End If

        Return

      End Subroutine confun
    End Module e04wb_a1w_fe_mod
    Program e04wb_a1w_fe

!     E04WB_A1W_F Example Main Program

!     .. Use Statements ..
      Use e04wb_a1w_fe_mod, Only: confun, nin, nout, objfun
      Use iso_c_binding, Only: c_ptr
      Use nagad_library, Only: e04ur_a1w_f, e04us_a1w_f, e04wb_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_w_rtype,         &
                               x10aa_a1w_f, x10ab_a1w_f, Assignment (=)
      Use nag_library, Only: nag_wp, x04caf
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Type (c_ptr)                     :: ad_handle
      Type (nagad_a1w_w_rtype)         :: objf
      Integer                          :: i, ifail, inform, iter, lb, lda,     &
                                          ldcj, ldfj, ldr, liwork, lwork, m,   &
                                          n, nclin, ncnln, sda
!     .. Local Arrays ..
      Type (nagad_a1w_w_rtype), Allocatable :: a(:,:), bl(:), bu(:), c(:),     &
                                          cjac(:,:), clamda(:), f(:),          &
                                          fjac(:,:), r(:,:), rwsav(:),         &
                                          work(:), x(:), y(:)
      Type (nagad_a1w_w_rtype)         :: ruser(44)
      Real (Kind=nag_wp), Allocatable  :: ar(:,:), blr(:), bur(:), dr(:),      &
                                          xr(:), yr(:)
      Real (Kind=nag_wp)               :: rr(44)
      Integer, Allocatable             :: istate(:), iwork(:), iwsav(:)
      Integer                          :: iuser(1)
      Logical, Allocatable             :: lwsav(:)
      Character (80)                   :: cwsav(1)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: max
!     .. Executable Statements ..
      Write (nout,*) 'E04WB_A1W_F Example Program Results'

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

!     Read the computational mode: 1 = algorithmic, 2 = symbolic
      Read (nin,*) m, n
      Read (nin,*) nclin, ncnln
      liwork = 3*n + nclin + 2*ncnln
      lda = max(1,nclin)

      If (nclin>0) Then
        sda = n
      Else
        sda = 1
      End If

      ldcj = max(1,ncnln)

      ldfj = m
      ldr = n

      If (ncnln==0 .And. nclin>0) Then
        lwork = 2*n**2 + 20*n + 11*nclin + m*(n+3)
      Else If (ncnln>0 .And. nclin>=0) Then
        lwork = 2*n**2 + n*nclin + 2*n*ncnln + 20*n + 11*nclin + 21*ncnln +    &
          m*(n+3)
      Else
        lwork = 20*n + m*(n+3)
      End If

      lb = n + nclin + ncnln
      Allocate (istate(lb),iwork(liwork),a(lda,sda),bl(lb),bu(lb),y(m),        &
        c(ncnln),cjac(ncnln,n),f(m),fjac(m,n),clamda(lb),r(ldr,n),x(n),        &
        work(lwork),lwsav(120),iwsav(610),rwsav(475),dr(m),yr(m),blr(lb),      &
        bur(lb),xr(n),ar(lda,sda))

      If (nclin>0) Then
        Read (nin,*)(ar(i,1:sda),i=1,nclin)
      End If
      a(1:nclin,1:sda) = ar(1:nclin,1:sda)

      Read (nin,*) yr(1:m)
      Read (nin,*) blr(1:lb)
      Read (nin,*) bur(1:lb)
      Read (nin,*) xr(1:n)
      y(1:m) = yr(1:m)
      bl(1:lb) = blr(1:lb)
      bu(1:lb) = bur(1:lb)
      x(1:n) = xr(1:n)

      rr(1:44) = (/8.0E0_nag_wp,8.0E0_nag_wp,10.0E0_nag_wp,10.0E0_nag_wp,      &
        10.0E0_nag_wp,10.0E0_nag_wp,12.0E0_nag_wp,12.0E0_nag_wp,12.0E0_nag_wp, &
        12.0E0_nag_wp,14.0E0_nag_wp,14.0E0_nag_wp,14.0E0_nag_wp,16.0E0_nag_wp, &
        16.0E0_nag_wp,16.0E0_nag_wp,18.0E0_nag_wp,18.0E0_nag_wp,20.0E0_nag_wp, &
        20.0E0_nag_wp,20.0E0_nag_wp,22.0E0_nag_wp,22.0E0_nag_wp,22.0E0_nag_wp, &
        24.0E0_nag_wp,24.0E0_nag_wp,24.0E0_nag_wp,26.0E0_nag_wp,26.0E0_nag_wp, &
        26.0E0_nag_wp,28.0E0_nag_wp,28.0E0_nag_wp,30.0E0_nag_wp,30.0E0_nag_wp, &
        30.0E0_nag_wp,32.0E0_nag_wp,32.0E0_nag_wp,34.0E0_nag_wp,36.0E0_nag_wp, &
        36.0E0_nag_wp,38.0E0_nag_wp,38.0E0_nag_wp,40.0E0_nag_wp,               &
        42.0E0_nag_wp/)

      ruser(1:44) = rr(1:44)

!     Create AD tape
      Call nagad_a1w_ir_create

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

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

!     Initialize sav arrays
      ifail = 0
      Call e04wb_a1w_f('E04USA',cwsav,1,lwsav,120,iwsav,610,rwsav,475,ifail)

!     Set option via string
      Call e04ur_a1w_f('Print Level = 1',lwsav,iwsav,rwsav,inform)

!     Solve the problem
      ifail = 0
      Call e04us_a1w_f(ad_handle,m,n,nclin,ncnln,lda,ldcj,ldfj,ldr,a,bl,bu,y,  &
        confun,objfun,iter,istate,c,cjac,f,fjac,clamda,objf,r,x,iwork,liwork,  &
        work,lwork,iuser,ruser,lwsav,iwsav,rwsav,ifail)

      Write (nout,99999) ' Optimal solution = ', objf%value
99999 Format (1X,A,F10.5)
      Write (nout,*)

      xr(1:n) = x(1:n)%value
      Call x04caf('General',' ',1,n,xr,1,' Solution point, x',ifail)

!     Primal results are printed by default

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

      Write (nout,*)
      Write (nout,*) ' Derivatives:'

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

!     Get derivatives
      dr(1:44) = nagad_a1w_get_derivative(ruser(1:44))

      Call x04caf('General',' ',1,44,dr,1,'     dx(1)/druser(1:44)',ifail)

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

    End Program e04wb_a1w_fe