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

NAG FL Interface Introduction
Example description
!   D01RLF Example Program Text
!   Mark 30.0 Release. NAG Copyright 2024.

    Module d01rlfe_mod

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

!     .. Use Statements ..
      Use nag_library, Only: nag_wp
!     .. Implicit None Statement ..
      Implicit None
!     .. Accessibility Statements ..
      Private
      Public                           :: f
!     .. Parameters ..
      Integer, Parameter, Public       :: nout = 6
    Contains
      Subroutine f(x,nx,fv,iflag,iuser,ruser,cpuser)

!       .. Use Statements ..
        Use, Intrinsic                 :: iso_c_binding, Only: c_ptr
!       .. Implicit None Statement ..
        Implicit None
!       .. Scalar Arguments ..
        Type (c_ptr), Intent (In)      :: cpuser
        Integer, Intent (Inout)        :: iflag
        Integer, Intent (In)           :: nx
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out) :: fv(nx)
        Real (Kind=nag_wp), Intent (Inout) :: ruser(*)
        Real (Kind=nag_wp), Intent (In) :: x(nx)
        Integer, Intent (Inout)        :: iuser(*)
!       .. Local Scalars ..
        Integer                        :: k
!       .. Intrinsic Procedures ..
        Intrinsic                      :: abs, any, sqrt
!       .. Executable Statements ..
        fv = abs(x-1.0E0_nag_wp/7.0E0_nag_wp)

        If (any(fv==0.0E0_nag_wp)) Then
!          A singular point will be hit.
!          Record offending abscissae and abort computation.
          iflag = 0
          Do k = 1, nx
            If (fv(k)==0.0E0_nag_wp) Then
              iflag = iflag + 1
              ruser(iflag) = x(k)
            End If
          End Do
!         store value of iflag in IUSER
          iuser(1) = iflag
!         signal abort by setting iflag<0
          iflag = -iflag

        Else
!         Safe to evaluate.
          fv(1:nx) = 1.0E0_nag_wp/sqrt(fv(1:nx))
        End If
        Return

      End Subroutine f
    End Module d01rlfe_mod
    Program d01rlfe

!     D01RLF Example Main Program

!     .. Use Statements ..
      Use d01rlfe_mod, Only: f, nout
      Use, Intrinsic                   :: iso_c_binding, Only: c_null_ptr,     &
                                          c_ptr
      Use nag_library, Only: d01rlf, nag_wp
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Type (c_ptr)                     :: cpuser
      Real (Kind=nag_wp)               :: a, abserr, b, epsabs, epsrel, result
      Integer                          :: ifail, liinfo, lrinfo, maxsub, npts
!     .. Local Arrays ..
      Real (Kind=nag_wp), Allocatable  :: points(:), rinfo(:), ruser(:)
      Integer, Allocatable             :: iinfo(:), iuser(:)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: max
!     .. Executable Statements ..
      Write (nout,*) 'D01RLF Example Program Results'

      epsabs = 0.0E0_nag_wp
      epsrel = 1.0E-04_nag_wp
      a = 0.0E0_nag_wp
      b = 1.0E0_nag_wp
      npts = 1
      maxsub = 20
      liinfo = 2*(max(maxsub,npts)) + npts + 4
      lrinfo = 4*(max(maxsub,npts)) + npts + 6

      Allocate (points(npts),rinfo(lrinfo),iinfo(liinfo),iuser(1),ruser(21))

      points(1) = 1.0E0_nag_wp/7.0E0_nag_wp
      iuser = 0
      ruser = 0.0E0_nag_wp
      cpuser = c_null_ptr

      ifail = -1
      Call d01rlf(f,a,b,npts,points,epsabs,epsrel,maxsub,result,abserr,rinfo,  &
        iinfo,iuser,ruser,cpuser,ifail)

      If (ifail>=0) Then
        Write (nout,*)
        Write (nout,99999) 'A       ', 'lower limit of integration', a
        Write (nout,99999) 'B       ', 'upper limit of integration', b
        Write (nout,99995) 'POINTS(1)', 'given break-point', points(1)
        Write (nout,99998) 'EPSABS  ', 'absolute accuracy requested', epsabs
        Write (nout,99998) 'EPSREL  ', 'relative accuracy requested', epsrel
        Write (nout,99996) 'MAXSUB  ', 'maximum number of subintervals',       &
          maxsub
        If (ifail<=5) Then
          Write (nout,*)
          Write (nout,99997) 'RESULT  ', 'approximation to the integral',      &
            result
          Write (nout,99998) 'ABSERR  ', 'estimate of the absolute error',     &
            abserr
          Write (nout,99996) 'IINFO(1)', 'number of subintervals used',        &
            iinfo(1)
        End If
      Else If (ifail==-1) Then
!       User requested exit.
        Write (nout,99994) ' Exit requested from F'
      End If

99999 Format (1X,A8,'  - ',A32,' = ',F9.4)
99998 Format (1X,A8,'  - ',A32,' = ',E9.2)
99997 Format (1X,A8,'  - ',A32,' = ',F9.5)
99996 Format (1X,A8,'  - ',A32,' = ',I4)
99995 Format (1X,A9,' - ',A32,' = ',F9.4)
99994 Format (1X,A30)

    End Program d01rlfe