program main !*****************************************************************************80 ! !! rbf_interp_1d_test() tests rbf_interp_1d(). ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 03 July 2018 ! ! Author: ! ! John Burkardt ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) integer nd external phi1 external phi2 external phi3 external phi4 external phi5 external phi6 external phi7 integer prob integer prob_num real ( kind = rk ) r0 real ( kind = rk ), allocatable :: xd(:) real ( kind = rk ) xmax real ( kind = rk ) xmin real ( kind = rk ), allocatable :: xy(:,:) call timestamp ( ) write ( *, '(a)' ) '' write ( *, '(a)' ) 'RBF_INTERP_1D_TEST():' write ( *, '(a)' ) ' FORTRAN90 version' write ( *, '(a)' ) ' Test RBF_INTERP_1D().' write ( *, '(a)' ) ' The R8LIB library is needed.' write ( *, '(a)' ) ' The test needs the TEST_INTERP library.' call p00_prob_num ( prob_num ) do prob = 1, prob_num ! ! Determine an appropriate value of R0, the spacing parameter. ! call p00_data_num ( prob, nd ) allocate ( xy(1:2,1:nd) ) call p00_data ( prob, 2, nd, xy ) allocate ( xd(1:nd) ) xd(1:nd) = xy(1,1:nd) xmax = maxval ( xd ) xmin = minval ( xd ) r0 = 3.0D+00 * ( xmax - xmin ) / real ( nd - 1, kind = rk ) deallocate ( xd ) deallocate ( xy ) call test01 ( prob, phi1, 'phi1', r0 ) call test01 ( prob, phi2, 'phi2', r0 ) call test01 ( prob, phi3, 'phi3', r0 ) call test01 ( prob, phi4, 'phi4', r0 ) call test01 ( prob, phi5, 'phi5', r0 ) call test01 ( prob, phi6, 'phi6', r0 ) call test01 ( prob, phi7, 'phi7', r0 ) end do ! ! Terminate. ! write ( *, '(a)' ) '' write ( *, '(a)' ) 'RBF_INTERP_1D_TEST:' write ( *, '(a)' ) ' Normal end of execution.' write ( *, '(a)' ) '' call timestamp ( ) stop end subroutine test01 ( prob, phi, phi_name, r0 ) !*****************************************************************************80 ! !! TEST01 tests RBF_INTERP_1D. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 18 September 2012 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer PROB, the index of the problem. ! ! Input, external PHI, the name of the radial basis function. ! ! Input, character ( len = * ) PHI_NAME, the name of the radial basis function. ! ! Input, real ( kind = rk ) R0, the scale factor. Typically, this might be ! a small multiple of the average distance between points. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) real ( kind = rk ) int_error real ( kind = rk ) ld real ( kind = rk ) li integer m integer nd integer ni external phi character ( len = * ) phi_name integer prob real ( kind = rk ) r0 real ( kind = rk ) r8vec_norm_affine real ( kind = rk ), allocatable :: w(:) real ( kind = rk ), allocatable :: xd(:) real ( kind = rk ), allocatable :: xi(:) real ( kind = rk ) xmax real ( kind = rk ) xmin real ( kind = rk ), allocatable :: xy(:,:) real ( kind = rk ), allocatable :: yd(:) real ( kind = rk ), allocatable :: yi(:) real ( kind = rk ) ymax real ( kind = rk ) ymin write ( *, '(a)' ) '' write ( *, '(a)' ) 'TEST01:' write ( *, '(a,i4)' ) ' Interpolate data from TEST_INTERP problem #', prob write ( *, '(a)' ) ' using radial basis function "' // trim ( phi_name ) // '".' write ( *, '(a,g14.6)' ) ' Scale factor R0 = ', r0 call p00_data_num ( prob, nd ) write ( *, '(a,i4)' ) ' Number of data points = ', nd allocate ( xy(1:2,1:nd) ) call p00_data ( prob, 2, nd, xy ) call r8mat_transpose_print ( 2, nd, xy, ' Data array:' ) allocate ( xd(1:nd) ) allocate ( yd(1:nd) ) xd(1:nd) = xy(1,1:nd) yd(1:nd) = xy(2,1:nd) m = 1 allocate ( w(1:nd) ) call rbf_weight ( m, nd, xd, r0, phi, yd, w ) ! ! #1: Does interpolant match function at interpolation points? ! ni = nd allocate ( xi(1:ni) ) allocate ( yi(1:ni) ) xi(1:nd) = xd(1:nd) call rbf_interp ( m, nd, xd, r0, phi, w, ni, xi, yi ) int_error = r8vec_norm_affine ( ni, yi, yd ) / real ( ni, kind = rk ) write ( *, '(a)' ) '' write ( *, '(a,g14.6)' ) & ' L2 interpolation error averaged per interpolant node = ', int_error deallocate ( xi ) deallocate ( yi ) ! ! #2: Compare estimated curve length to piecewise linear (minimal) curve length. ! Assume data is sorted, and normalize X and Y dimensions by (XMAX-XMIN) and ! (YMAX-YMIN). ! xmin = minval ( xd(1:nd) ) xmax = maxval ( xd(1:nd) ) ymin = minval ( yd(1:nd) ) ymax = maxval ( yd(1:nd) ) ni = 501 allocate ( xi(1:ni) ) allocate ( yi(1:ni) ) call r8vec_linspace ( ni, xmin, xmax, xi ) call rbf_interp ( m, nd, xd, r0, phi, w, ni, xi, yi ) ld = sum ( sqrt ( ( ( xd(2:nd) - xd(1:nd-1) ) / ( xmax - xmin ) ) ** 2 & + ( ( yd(2:nd) - yd(1:nd-1) ) / ( ymax - ymin ) ) ** 2 ) ) li = sum ( sqrt ( ( ( xi(2:ni) - xi(1:ni-1) ) / ( xmax - xmin ) ) ** 2 & + ( ( yi(2:ni) - yi(1:ni-1) ) / ( ymax - ymin ) ) ** 2 ) ) write ( *, '(a)' ) '' write ( *, '(a,g14.6)' ) ' Normalized length of piecewise linear interpolant = ', ld write ( *, '(a,g14.6)' ) ' Normalized length of polynomial interpolant = ', li deallocate ( w ) deallocate ( xd ) deallocate ( xi ) deallocate ( xy ) deallocate ( yd ) deallocate ( yi ) return end subroutine timestamp ( ) !*****************************************************************************80 ! !! timestamp() prints the current YMDHMS date as a time stamp. ! ! Example: ! ! 31 May 2001 9:45:54.872 AM ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 15 August 2021 ! ! Author: ! ! John Burkardt ! implicit none character ( len = 8 ) ampm integer d integer h integer m integer mm character ( len = 9 ), parameter, dimension(12) :: month = (/ & 'January ', 'February ', 'March ', 'April ', & 'May ', 'June ', 'July ', 'August ', & 'September', 'October ', 'November ', 'December ' /) integer n integer s integer values(8) integer y call date_and_time ( values = values ) y = values(1) m = values(2) d = values(3) h = values(5) n = values(6) s = values(7) mm = values(8) if ( h < 12 ) then ampm = 'AM' else if ( h == 12 ) then if ( n == 0 .and. s == 0 ) then ampm = 'Noon' else ampm = 'PM' end if else h = h - 12 if ( h < 12 ) then ampm = 'PM' else if ( h == 12 ) then if ( n == 0 .and. s == 0 ) then ampm = 'Midnight' else ampm = 'AM' end if end if end if write ( *, '(i2.2,1x,a,1x,i4,2x,i2,a1,i2.2,a1,i2.2,a1,i3.3,1x,a)' ) & d, trim ( month(m) ), y, h, ':', n, ':', s, '.', mm, trim ( ampm ) return end