program main c*********************************************************************72 c cc vandermonde_interp_1d_test() tests vandermonde_interp_1d(). c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 21 September 2012 c c Author: c c John Burkardt c implicit none integer prob integer prob_num call timestamp ( ) write ( *, '(a)' ) '' write ( *, '(a)' ) 'VANDERMONDE_INTERP_1D_test():' write ( *, '(a)' ) ' Fortran77 version' write ( *, '(a)' ) ' Test the VANDERMONDE_INTERP_1D library.' write ( *, '(a)' ) ' The R8LIB library is needed.' write ( *, '(a)' ) ' The QR_SOLVE library is needed.' write ( *, '(a)' ) ' This test needs the TEST_INTERP library.' write ( *, '(a)' ) ' This test needs the CONDITION library.' call p00_prob_num ( prob_num ) do prob = 1, prob_num call test01 ( prob ) end do do prob = 1, prob_num call test02 ( prob ) end do c c Terminate. c write ( *, '(a)' ) '' write ( *, '(a)' ) 'VANDERMONDE_INTERP_1D_test():' write ( *, '(a)' ) ' Normal end of execution.' write ( *, '(a)' ) '' call timestamp ( ) return end subroutine test01 ( prob ) c*********************************************************************72 c cc TEST01 tests VANDERMONDE_INTERP_1D_MATRIX. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 07 October 2012 c c Author: c c John Burkardt c implicit none integer nd_max parameter ( nd_max = 49 ) integer ni_max parameter ( ni_max = 501 ) double precision ad(nd_max,nd_max) double precision cd(nd_max) double precision condition logical debug parameter ( debug = .false. ) integer i double precision int_error double precision ld double precision li integer m integer nd integer ni integer prob double precision r8vec_max double precision r8vec_min double precision r8vec_norm_affine double precision xd(nd_max) double precision xi(ni_max) double precision xmax double precision xmin double precision xy(2,nd_max) double precision yd(nd_max) double precision yi(ni_max) double precision ymax double precision ymin write ( *, '(a)' ) '' write ( *, '(a)' ) 'TEST01:' write ( *, '(a,i2)' ) & ' Interpolate data from TEST_INTERP problem #', prob call p00_data_num ( prob, nd ) write ( *, '(a,i2)' ) ' Number of data points = ', nd call p00_data ( prob, 2, nd, xy ) if ( debug ) then call r8mat_transpose_print ( 2, nd, xy, ' Data array:' ) end if do i = 1, nd xd(i) = xy(1,i) yd(i) = xy(2,i) end do c c Choose the degree of the polynomial to be ND - 1. c m = nd - 1 c c Compute Vandermonde matrix and get condition number. c call vandermonde_interp_1d_matrix ( nd, xd, ad ) call condition_hager ( nd, ad, condition ) write ( *, '(a)' ) '' write ( *, '(a,g14.6)' ) & ' Condition of Vandermonde matrix is ', condition c c Solve linear system. c call qr_solve ( nd, nd, ad, yd, cd ) c c #1: Does interpolant match function at interpolation points? c ni = nd do i = 1, ni xi(i) = xd(i) end do call vandermonde_value_1d ( nd, cd, ni, xi, yi ) int_error = r8vec_norm_affine ( ni, yi, yd ) / dble ( ni ) write ( *, '(a)' ) '' write ( *, '(a,g14.6)' ) & ' L2 interpolation error averaged per interpolant node = ', & int_error c c #2: Compare estimated curve length to piecewise linear (minimal) curve length. c Assume data is sorted, and normalize X and Y dimensions by (XMAX-XMIN) and c (YMAX-YMIN). c xmin = r8vec_min ( nd, xd ) xmax = r8vec_max ( nd, xd ) ymin = r8vec_min ( nd, yd ) ymax = r8vec_max ( nd, yd ) ni = 501 call r8vec_linspace ( ni, xmin, xmax, xi ) call vandermonde_value_1d ( nd, cd, ni, xi, yi ) ld = 0.0D+00 do i = 1, nd - 1 ld = ld + sqrt & ( ( ( xd(i+1) - xd(i) ) / ( xmax - xmin ) ) ** 2 & + ( ( yd(i+1) - yd(i) ) / ( ymax - ymin ) ) ** 2 ) end do li = 0.0D+00 do i = 1, ni - 1 li = li + sqrt & ( ( ( xi(i+1) - xi(i) ) / ( xmax - xmin ) ) ** 2 & + ( ( yi(i+1) - yi(i) ) / ( ymax - ymin ) ) ** 2 ) end do write ( *, '(a)' ) '' write ( *, '(a,g14.6)' ) & ' Normalized length of piecewise linear interpolant = ', ld write ( *, '(a,g14.6)' ) & ' Normalized length of polynomial interpolant = ', li return end subroutine test02 ( prob ) c*********************************************************************72 c cc TEST02 tests VANDERMONDE_INTERP_1D_MATRIX. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 01 June 2013 c c Author: c c John Burkardt c c Parameters: c c Input, integer PROB, the problem index. c implicit none integer nd_max parameter ( nd_max = 49 ) integer ni_max parameter ( ni_max = 501 ) double precision ad(nd_max,nd_max) double precision cd(nd_max) character * ( 255 ) command_filename integer command_unit character * ( 255 ) data_filename integer data_unit integer i character * ( 255 ) interp_filename integer interp_unit integer j integer nd integer ni character * ( 255 ) output_filename integer prob double precision r8vec_max double precision r8vec_min double precision xd(nd_max) double precision xi(ni_max) double precision xmax double precision xmin double precision xy(2,nd_max) double precision yd(nd_max) double precision yi(ni_max) write ( *, '(a)' ) '' write ( *, '(a)' ) 'TEST02:' write ( *, '(a)' ) & ' VANDERMONDE_INTERP_1D_MATRIX sets the Vandermonde linear' write ( *, '(a)' ) ' system for the interpolating polynomial.' write ( *, '(a,i2)' ) & ' Interpolate data from TEST_INTERP problem #', prob call p00_data_num ( prob, nd ) write ( *, '(a,i4)' ) ' Number of data points = ', nd call p00_data ( prob, 2, nd, xy ) call r8mat_transpose_print ( 2, nd, xy, ' Data array:' ) do i = 1, nd xd(i) = xy(1,i) yd(i) = xy(2,i) end do c c Set up the Vandermonde matrix AD. c call vandermonde_interp_1d_matrix ( nd, xd, ad ) c c Solve the linear system for the polynomial coefficients CD. c call qr_solve ( nd, nd, ad, yd, cd ) c c Create data file. c write ( data_filename, '(a,i2.2,a)' ) 'data', prob, '.txt' call get_unit ( data_unit ) open ( unit = data_unit, file = data_filename, & status = 'replace' ) do j = 1, nd write ( data_unit, '(2x,g14.6,2x,g14.6)' ) xd(j), yd(j) end do close ( unit = data_unit ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) & ' Created graphics data file "' & // trim ( data_filename ) // '".' c c Create interp file. c ni = 501 xmin = r8vec_min ( nd, xd ) xmax = r8vec_max ( nd, xd ) call r8vec_linspace ( ni, xmin, xmax, xi ) call vandermonde_value_1d ( nd, cd, ni, xi, yi ) write ( interp_filename, '(a,i2.2,a)' ) 'interp', prob, '.txt' call get_unit ( interp_unit ) open ( unit = interp_unit, file = interp_filename, & status = 'replace' ) do j = 1, ni write ( interp_unit, '(2x,g14.6,2x,g14.6)' ) xi(j), yi(j) end do close ( unit = interp_unit ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) & ' Created graphics interp file "' & // trim ( interp_filename ) // '".' c c Plot the data and the interpolant. c write ( command_filename, '(a,i2.2,a)' ) 'commands', prob, '.txt' call get_unit ( command_unit ) open ( unit = command_unit, file = command_filename, & status = 'replace' ) write ( output_filename, '(a,i2.2,a)' ) 'plot', prob, '.png' write ( command_unit, '(a)' ) '# ' // trim ( command_filename ) write ( command_unit, '(a)' ) '#' write ( command_unit, '(a)' ) '# Usage:' write ( command_unit, '(a)' ) & '# gnuplot < ' // trim ( command_filename ) write ( command_unit, '(a)' ) '#' write ( command_unit, '(a)' ) 'set term png' write ( command_unit, '(a)' ) & 'set output "' // trim ( output_filename ) // '"' write ( command_unit, '(a)' ) 'set xlabel "<---X--->"' write ( command_unit, '(a)' ) 'set ylabel "<---Y--->"' write ( command_unit, '(a)' ) & 'set title "Data versus Vandermonde Polynomial Interpolant"' write ( command_unit, '(a)' ) 'set grid' write ( command_unit, '(a)' ) 'set style data lines' write ( command_unit, '(a)' ) & 'plot "' // trim ( data_filename ) // & '" using 1:2 with points pt 7 ps 2 lc rgb "blue",\' write ( command_unit, '(a)' ) & ' "' // trim ( interp_filename ) // & '" using 1:2 lw 3 linecolor rgb "red"' close ( unit = command_unit ) write ( *, '(a)' ) & ' Created graphics command file "' & // trim ( command_filename ) // '".' return end subroutine get_unit ( iunit ) c*********************************************************************72 c cc get_unit() returns a free Fortran unit number. c c Discussion: c c A "free" Fortran unit number is a value between 1 and 99 which c is not currently associated with an I/O device. A free Fortran unit c number is needed in order to open a file with the OPEN command. c c If IUNIT = 0, then no free Fortran unit could be found, although c all 99 units were checked (except for units 5, 6 and 9, which c are commonly reserved for console I/O). c c Otherwise, IUNIT is a value between 1 and 99, representing a c free Fortran unit. Note that GET_UNIT assumes that units 5 and 6 c are special, and will never return those values. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 02 September 2013 c c Author: c c John Burkardt c c Output: c c integer IUNIT, the free unit number. c implicit none integer i integer iunit logical value iunit = 0 do i = 1, 99 if ( i .ne. 5 .and. i .ne. 6 .and. i .ne. 9 ) then inquire ( unit = i, opened = value, err = 10 ) if ( .not. value ) then iunit = i return end if end if 10 continue end do return end subroutine r8mat_transpose_print ( m, n, a, title ) c*********************************************************************72 c cc r8mat_transpose_print() prints an R8MAT, transposed. c c Discussion: c c An R8MAT is an array of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 28 April 2008 c c Author: c c John Burkardt c c Parameters: c c Input, integer M, N, the number of rows and columns. c c Input, double precision A(M,N), an M by N matrix to be printed. c c Input, character*(*) TITLE, a title. c implicit none integer m integer n double precision a(m,n) character*(*) title call r8mat_transpose_print_some ( m, n, a, 1, 1, m, n, title ) return end subroutine r8mat_transpose_print_some ( m, n, a, ilo, jlo, ihi, & jhi, title ) c*********************************************************************72 c cc r8mat_transpose_print_some() prints some of an R8MAT transposed. c c Discussion: c c An R8MAT is an array of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 28 April 2008 c c Author: c c John Burkardt c c Parameters: c c Input, integer M, N, the number of rows and columns. c c Input, double precision A(M,N), an M by N matrix to be printed. c c Input, integer ILO, JLO, the first row and column to print. c c Input, integer IHI, JHI, the last row and column to print. c c Input, character * ( * ) TITLE, a title. c implicit none integer incx parameter ( incx = 5 ) integer m integer n double precision a(m,n) character * ( 14 ) ctemp(incx) integer i integer i2 integer i2hi integer i2lo integer ihi integer ilo integer inc integer j integer j2hi integer j2lo integer jhi integer jlo character * ( * ) title write ( *, '(a)' ) ' ' write ( *, '(a)' ) trim ( title ) if ( m .le. 0 .or. n .le. 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' (None)' return end if do i2lo = max ( ilo, 1 ), min ( ihi, m ), incx i2hi = i2lo + incx - 1 i2hi = min ( i2hi, m ) i2hi = min ( i2hi, ihi ) inc = i2hi + 1 - i2lo write ( *, '(a)' ) ' ' do i = i2lo, i2hi i2 = i + 1 - i2lo write ( ctemp(i2), '(i8,6x)') i end do write ( *, '('' Row'',5a14)' ) ctemp(1:inc) write ( *, '(a)' ) ' Col' j2lo = max ( jlo, 1 ) j2hi = min ( jhi, n ) do j = j2lo, j2hi do i2 = 1, inc i = i2lo - 1 + i2 write ( ctemp(i2), '(g14.6)' ) a(i,j) end do write ( *, '(2x,i8,a,5a14)' ) j, ':', ( ctemp(i), i = 1, inc ) end do end do return end subroutine r8vec_linspace ( n, a, b, x ) c*********************************************************************72 c cc r8vec_linspace() creates a vector of linearly spaced values. c c Discussion: c c An R8VEC is a vector of R8's. c c 4 points evenly spaced between 0 and 12 will yield 0, 4, 8, 12. c c In other words, the interval is divided into N-1 even subintervals, c and the endpoints of intervals are used as the points. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 14 March 2011 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the number of entries in the vector. c c Input, double precision A, B, the first and last entries. c c Output, double precision X(N), a vector of linearly spaced data. c implicit none integer n double precision a double precision b integer i double precision x(n) if ( n .eq. 1 ) then x(1) = ( a + b ) / 2.0D+00 else do i = 1, n x(i) = ( dble ( n - i ) * a & + dble ( i - 1 ) * b ) & / dble ( n - 1 ) end do end if return end function r8vec_max ( n, a ) c*********************************************************************72 c cc r8vec_max() returns the maximum value in an R8VEC. c c Discussion: c c An R8VEC is a vector of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 12 May 2014 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the number of entries in the array. c c Input, double precision A(N), the array. c c Output, double precision R8VEC_MAX, the value of the largest entry. c implicit none integer n double precision a(n) integer i double precision r8_huge parameter ( r8_huge = 1.79769313486231571D+308 ) double precision r8vec_max double precision value value = - r8_huge do i = 1, n value = max ( value, a(i) ) end do r8vec_max = value return end function r8vec_min ( n, a ) c*********************************************************************72 c cc r8vec_min() returns the minimum value in an R8VEC. c c Discussion: c c An R8VEC is a vector of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 18 January 2015 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the number of entries in the array. c c Input, double precision A(N), the array. c c Output, double precision R8VEC_MIN, the value of the smallest entry. c implicit none integer n double precision a(n) integer i double precision r8_huge parameter ( r8_huge = 1.79769313486231571D+308 ) double precision r8vec_min double precision value value = r8_huge do i = 1, n value = min ( value, a(i) ) end do r8vec_min = value return end function r8vec_norm_affine ( n, v0, v1 ) c*********************************************************************72 c cc r8vec_norm_affine() returns the affine norm of an R8VEC. c c Discussion: c c An R8VEC is a vector of R8's. c c The affine vector L2 norm is defined as: c c R8VEC_NORM_AFFINE(V0,V1) c = sqrt ( sum ( 1 <= I <= N ) ( V1(I) - V0(I) )^2 ) c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 27 October 2010 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the order of the vectors. c c Input, double precision V0(N), the base vector. c c Input, double precision V1(N), the vector whose affine norm is desired. c c Output, double precision R8VEC_NORM_AFFINE, the L2 norm of V1-V0. c implicit none integer n integer i double precision r8vec_norm_affine double precision v0(n) double precision v1(n) r8vec_norm_affine = 0.0D+00 do i = 1, n r8vec_norm_affine = r8vec_norm_affine & + ( v0(i) - v1(i) ) ** 2 end do r8vec_norm_affine = sqrt ( r8vec_norm_affine ) return end subroutine timestamp ( ) c*********************************************************************72 c cc timestamp() prints the YMDHMS date as a timestamp. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 12 June 2014 c c Author: c c John Burkardt c implicit none character * ( 8 ) ampm integer d character * ( 8 ) date integer h integer m integer mm character * ( 9 ) month(12) integer n integer s character * ( 10 ) time integer y save month data month / & 'January ', 'February ', 'March ', 'April ', & 'May ', 'June ', 'July ', 'August ', & 'September', 'October ', 'November ', 'December ' / call date_and_time ( date, time ) read ( date, '(i4,i2,i2)' ) y, m, d read ( time, '(i2,i2,i2,1x,i3)' ) h, n, s, mm if ( h .lt. 12 ) then ampm = 'AM' else if ( h .eq. 12 ) then if ( n .eq. 0 .and. s .eq. 0 ) then ampm = 'Noon' else ampm = 'PM' end if else h = h - 12 if ( h .lt. 12 ) then ampm = 'PM' else if ( h .eq. 12 ) then if ( n .eq. 0 .and. s .eq. 0 ) then ampm = 'Midnight' else ampm = 'AM' end if end if end if write ( *, & '(i2,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