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

NAG AD Library Introduction
Example description
    Program g02da_t1w_fe

!     G02DA_T1W_F Example Program Text
!     Mark 30.1 Release. NAG Copyright 2024.

!     .. Use Statements ..
      Use iso_c_binding, Only: c_ptr
      Use nagad_library, Only: g02da_t1w_f, nagad_t1w_get, nagad_t1w_set,      &
                               nagad_t1w_w_rtype, x10aa_t1w_f, x10ab_t1w_f
      Use nag_library, Only: nag_wp, x04caf
!     .. Implicit None Statement ..
      Implicit None
!     .. Parameters ..
      Real (Kind=nag_wp), Parameter    :: zero = 0.0_nag_wp
      Integer, Parameter               :: nin = 5, nout = 6
!     .. Local Scalars ..
      Type (c_ptr)                     :: ad_handle
      Type (nagad_t1w_w_rtype)         :: rss, tol
      Real (Kind=nag_wp)               :: sav, tlm
      Integer                          :: i, idf, ifail, ip, irank, j, ldq,    &
                                          ldx, lwt, m, n
      Logical                          :: svd
      Character (1)                    :: mean, weight
!     .. Local Arrays ..
      Type (nagad_t1w_w_rtype), Allocatable :: b(:), cov(:), h(:), p(:),       &
                                          q(:,:), res(:), se(:), wk(:), wt(:), &
                                          x(:,:), y(:)
      Real (Kind=nag_wp), Allocatable  :: dbdy(:,:), dy(:)
      Integer, Allocatable             :: isx(:)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: count
!     .. Executable Statements ..
      Write (nout,*) 'G02DA_T1W_F Example Program Results'
      Flush (nout)

!     Skip heading in data file
      Read (nin,*)
      Read (nin,*) n, m, weight, mean

      If (weight=='W' .Or. weight=='w') Then
        lwt = n
      Else
        lwt = 0
      End If
      ldx = n
      Allocate (x(ldx,m),y(n),wt(lwt),isx(m),dy(n))

!     Read in data
      x(1:n,1:m)%tangent = zero
      y(1:n)%tangent = zero
      wt(1:lwt)%tangent = zero
      If (lwt>0) Then
        Read (nin,*)(x(i,1:m)%value,y(i)%value,wt(i)%value,i=1,n)
      Else
        Read (nin,*)(x(i,1:m)%value,y(i)%value,i=1,n)
      End If

!     Read in variable inclusion flags
      Read (nin,*) isx(1:m)

!     Calculate IP
      ip = count(isx(1:m)>0)
      If (mean=='M' .Or. mean=='m') Then
        ip = ip + 1
      End If

      ldq = n

      Allocate (b(ip),cov((ip*ip+ip)/2),h(n),p(ip*(ip+                         &
        2)),q(ldq,ip+1),res(n),se(ip),wk(ip*ip+5*(ip-1)),dbdy(n,ip))

!     Use suggested value for tolerance
      tol%value = 0.000001E0_nag_wp
      tol%tangent = zero

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

!     First call for primal results
      ifail = -1
      Call g02da_t1w_f(ad_handle,mean,weight,n,x,ldx,m,isx,ip,y,wt,rss,idf,b,  &
        se,cov,res,h,q,ldq,svd,irank,p,tol,wk,ifail)
      If (ifail/=0) Then
        If (ifail/=5) Then
          Go To 100
        End If
      End If
!     Display primal results on first call
      If (svd) Then
        Write (nout,99999) 'Model not of full rank, rank = ', irank
        Write (nout,*)
      End If
      Write (nout,99998) 'Residual sum of squares = ', rss%value
      Write (nout,99999) 'Degrees of freedom      = ', idf
      Write (nout,*)
      Write (nout,*) 'Variable   Parameter estimate   Standard error'
      Write (nout,*)
      If (ifail==0) Then
        Write (nout,99997)(j,b(j)%value,se(j)%value,j=1,ip)
      Else
        Write (nout,99996)(j,b(j)%value,j=1,ip)
      End If

!     Repeat Tangent linear call for each input y(1:n)
      Do i = 1, n
        sav = y(i)%value
        Call nagad_t1w_set(y(i),1.0_nag_wp,1)

!       Fit general linear regression model
        ifail = -1
        Call g02da_t1w_f(ad_handle,mean,weight,n,x,ldx,m,isx,ip,y,wt,rss,idf,  &
          b,se,cov,res,h,q,ldq,svd,irank,p,tol,wk,ifail)

!       Get derivative information w.r.t. y(i)
        tlm = 0.0_nag_wp
        Call nagad_t1w_get(rss,tlm,1)
        dy(i) = tlm
        Do j = 1, ip
          tlm = 0.0_nag_wp
          Call nagad_t1w_get(b(j),tlm,1)
          dbdy(i,j) = tlm
        End Do

!       Reset input y(i)
        Call nagad_t1w_set(y(i),0.0_nag_wp,1)
        y(i)%value = sav
      End Do

!     Display derivative results
      Write (nout,*)
      Write (nout,*) ' Derivatives calculated: First order tangent linear'
      Write (nout,*) ' Computational mode    : algorithmic'

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

      Write (nout,*) '    i        d(rss)/dy(i) '
      Do i = 1, n
        Write (nout,99995) i, dy(i)
      End Do

      Write (nout,*)
      ifail = 0
      Call x04caf('General',' ',n,ip,dbdy,n,'    db/dy',ifail)

100   Continue

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


99999 Format (1X,A,I4)
99998 Format (1X,A,E12.4)
99997 Format (1X,I6,2E20.4)
99996 Format (1X,I6,E20.4)
99995 Format (1X,I5,6X,F9.4)
    End Program g02da_t1w_fe