Program f11jb_a1w_fe
! F11JB_A1W_F Example Program Text
! Mark 30.2 Release. NAG Copyright 2024.
! .. Use Statements ..
Use iso_c_binding, Only: c_ptr
Use nagad_library, Only: f11ja_a1w_f, f11jb_a1w_f, f11za_a1w_f, &
f11zb_a1w_f, nagad_a1w_get_derivative, &
nagad_a1w_inc_derivative, &
nagad_a1w_ir_interpret_adjoint_sparse, &
nagad_a1w_ir_register_variable, &
nagad_a1w_ir_remove, nagad_a1w_ir_zero_adjoints &
, nagad_a1w_w_rtype, x10aa_a1w_f, x10ab_a1w_f, &
x10za_a1w_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
Type (nagad_a1w_w_rtype) :: dscale, dtol
Integer :: i, ifail, j, la, lfill, liwork, n, &
nnz, nnzc, npivm
Character (1) :: check, mic, pstrat
! .. Local Arrays ..
Type (nagad_a1w_w_rtype), Allocatable :: a(:), x(:), y(:)
Real (Kind=nag_wp), Allocatable :: ar(:), dxdy(:,:), yr(:)
Integer, Allocatable :: icol(:), ipiv(:), irow(:), istr(:), &
iwork(:), perm_fwd(:), perm_inv(:)
! .. Executable Statements ..
Write (nout,*) 'F11JB_A1W_F Example Program Results'
! Skip heading in data file
Read (nin,*)
! Read order of matrix and number of nonzero entries
Read (nin,*) n
Read (nin,*) nnz
la = 3*nnz
liwork = 2*la + 7*n + 1
Allocate (a(la),x(n),y(n),icol(la),ipiv(n),irow(la),istr(n+1), &
iwork(liwork),perm_fwd(n),perm_inv(n))
Allocate (ar(la),yr(n),dxdy(n,n))
! Read the matrix A
Do i = 1, nnz
Read (nin,*) ar(i), irow(i), icol(i)
End Do
a(1:nnz) = ar(1:nnz)
! Read the vector y
Read (nin,*) yr(1:n)
y = yr
! Create AD tape
Call x10za_a1w_f
! Create AD configuration data object
ifail = 0
Call x10aa_a1w_f(ad_handle,ifail)
! Calculate Cholesky factorization
lfill = -1
dtol = 0.0E0_nag_wp
mic = 'N'
dscale = 0.0E0_nag_wp
pstrat = 'M'
! Register variables to differentiate w.r.t.
Call nagad_a1w_ir_register_variable(y)
iwork = 0
! Compute reverse Cuthill-McKee permutation for bandwidth reduction
Call do_rcm(ad_handle,irow,icol,a,y,istr,perm_fwd,perm_inv,iwork)
ifail = 0
Call f11ja_a1w_f(ad_handle,n,nnz,a,la,irow,icol,lfill,dtol,mic,dscale, &
pstrat,ipiv,istr,nnzc,npivm,iwork,liwork,ifail)
! Check the output value of NPIVM
If (npivm/=0) Then
Write (nout,99998) 'Factorization is not complete', npivm
Else
! Solve P L D L^T P^T x = y
x = 0.0_nag_wp
check = 'C'
ifail = 0
Call f11jb_a1w_f(ad_handle,n,a,la,irow,icol,ipiv,istr,check,y,x,ifail)
! Output results
Write (nout,*) ' Solution of linear system with Reverse Cuthill-McKee'
Write (nout,99999)(x(perm_inv(i))%value,i=1,n)
Write (nout,*)
Write (nout,*) ' Derivatives calculated: First order adjoints'
Write (nout,*) ' Computational mode : algorithmic'
Write (nout,*)
Write (nout,*) ' Derivatives of solution X w.r.t. RHS Y (A inverse)'
Write (nout,*)
! Setup evaluation of derivatives via adjoints
Do i = 1, n
Call nagad_a1w_ir_zero_adjoints
Call nagad_a1w_inc_derivative(x(perm_inv(i)),1.0_nag_wp)
ifail = 0
Call nagad_a1w_ir_interpret_adjoint_sparse(ifail)
! Get derivatives
Do j = 1, n
dxdy(i,j) = nagad_a1w_get_derivative(y(perm_inv(j)))
End Do
End Do
Call x04caf('General',' ',n,n,dxdy,n,' dx(i)/dy(j)',ifail)
End If
! Remove computational data object and tape
Call x10ab_a1w_f(ad_handle,ifail)
Call nagad_a1w_ir_remove
99999 Format (1X,E16.4)
99998 Format (1X,A,I20)
Contains
Subroutine do_rcm(ad_handle,irow,icol,a,y,istr,perm_fwd,perm_inv,iwork)
! .. Use Statements ..
Use iso_c_binding, Only: c_ptr
Use nag_library, Only: f11yef
! .. Parameters ..
Logical, Parameter :: lopts(5) = (/.False.,.False.,.True., &
.True.,.True./)
! .. Scalar Arguments ..
Type (c_ptr) :: ad_handle
! .. Array Arguments ..
Type (nagad_a1w_w_rtype), Intent (Inout) :: a(la), y(n)
Integer, Intent (Inout) :: icol(la), irow(la), istr(n+1), &
iwork(*)
Integer, Intent (Out) :: perm_fwd(n), perm_inv(n)
! .. Local Scalars ..
Integer :: i, ifail, j, nnz_cs, nnz_scs
! .. Local Arrays ..
Type (nagad_a1w_w_rtype), Allocatable :: rwork(:)
Integer :: info(4), mask(1)
! .. Intrinsic Procedures ..
Intrinsic :: size
! .. Executable Statements ..
! SCS to CS, must add the upper triangle entries.
j = nnz + 1
Do i = 1, nnz
If (irow(i)>icol(i)) Then
! strictly lower triangle, add the transposed
a(j) = a(i)
irow(j) = icol(i)
icol(j) = irow(i)
j = j + 1
End If
End Do
nnz_cs = j - 1
! Reorder, CS to CCS, icolzp in istr
ifail = 0
Call f11za_a1w_f(ad_handle,n,nnz_cs,a,icol,irow,'F','F',istr,iwork, &
ifail)
! Calculate reverse Cuthill-McKee
ifail = 0
Call f11yef(n,nnz_cs,istr,irow,lopts,mask,perm_fwd,info,ifail)
! compute inverse perm, in perm_inv(1:n)
Do i = 1, n
perm_inv(perm_fwd(i)) = i
End Do
! Apply permutation on column/row indices
icol(1:nnz_cs) = perm_inv(icol(1:nnz_cs))
irow(1:nnz_cs) = perm_inv(irow(1:nnz_cs))
! restrict to lower triangle, SCS format
! copying entries upwards
j = 1
Do i = 1, nnz_cs
If (irow(i)>=icol(i)) Then
! non-upper triangle, bubble up
a(j) = a(i)
icol(j) = icol(i)
irow(j) = irow(i)
j = j + 1
End If
End Do
nnz_scs = j - 1
! sort
ifail = 0
Call f11zb_a1w_f(ad_handle,n,nnz_scs,a,irow,icol,'S','K',istr,iwork, &
ifail)
! permute rhs vector
Allocate (rwork(size(perm_fwd)))
rwork(:) = y(perm_fwd(:))
y(:) = rwork(:)
Deallocate (rwork)
End Subroutine do_rcm
End Program f11jb_a1w_fe