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

NAG AD Library Introduction
Example description
    Program f07ad_t1w_fe

!     F07AD_T1W_F Example Program Text
!     Mark 30.0 Release. NAG Copyright 2024.

!     .. Use Statements ..
      Use iso_c_binding, Only: c_ptr
      Use nagad_library, Only: f07ad_t1w_f, f07aj_t1w_f,                       &
                               nagad_t1w_set_derivative, nagad_t1w_w_rtype,    &
                               x10aa_t1w_f, x10ab_t1w_f, Assignment (=)
      Use nag_library, Only: nag_wp, x04caf
!     .. Implicit None Statement ..
      Implicit None
!     .. Parameters ..
      Integer, Parameter               :: nin = 5, nout = 6
!     .. Local Scalars ..
      Type (c_ptr)                     :: ad_handle
      Integer                          :: i, ifail, j, lwork, n
!     .. Local Arrays ..
      Type (nagad_t1w_w_rtype), Allocatable :: a(:,:), a_in(:,:), work(:)
      Real (Kind=nag_wp), Allocatable  :: ar(:,:), da(:,:)
      Integer, Allocatable             :: ipiv(:)
!     .. Executable Statements ..
      Write (nout,*) 'F07AD_T1W_F Example Program Results'
!     Skip heading in data file
      Read (nin,*)
      Read (nin,*) n

      lwork = 64*n
      Allocate (a(n,n),a_in(n,n),ipiv(n))
      Allocate (ar(n,n),da(n,n),work(lwork))

!     Read A from data file
      Read (nin,*)(ar(i,1:n),i=1,n)
      a_in = ar

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

      Do i = 1, n
        Call nagad_t1w_set_derivative(a_in(i,i),1.0_nag_wp)
        a = a_in

!       Factorize A
        ifail = 0
        Call f07ad_t1w_f(ad_handle,n,n,a,n,ipiv,ifail)

!       Compute inverse
        Call f07aj_t1w_f(ad_handle,n,a,n,ipiv,work,lwork,ifail)

        If (i==1) Then
!         Print solution
          ar(1:n,1:n) = a(1:n,1:n)%value

          Write (nout,*)
          Flush (nout)

!         ifail: behaviour on error exit
!             =0 for hard exit, =1 for quiet-soft, =-1 for noisy-soft
          ifail = 0
          Call x04caf('General',' ',n,n,ar,n,'Inverse',ifail)
        End If
        Do j = 1, n
          da(j,i) = a(j,j)%tangent
        End Do
        a_in(i,i)%tangent = 0.0_nag_wp
      End Do

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

      Write (nout,*)
      Write (nout,*) ' Derivatives of inverse diagonal w.r.t. diagonal of A'
      Write (nout,*)

      Call x04caf('General',' ',n,n,da,n,'       d(inv(i,i))/da(j,j)',ifail)

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

    End Program f07ad_t1w_fe