program main !*****************************************************************************80 ! !! MAIN is the main program for JACOBI_RULE. ! ! Discussion: ! ! This program computes a standard Gauss-Jacobi quadrature rule ! and writes it to a file. ! ! The user specifies: ! * the ORDER (number of points) in the rule; ! * ALPHA and BETA parameters; ! * the root name of the output files. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 16 January 2008 ! ! Author: ! ! John Burkardt ! implicit none real ( kind = 8 ) alpha integer ( kind = 4 ) arg_num real ( kind = 8 ) beta integer ( kind = 4 ) iarg integer ( kind = 4 ) iargc integer ( kind = 4 ) ierror integer ( kind = 4 ) last integer ( kind = 4 ) order character ( len = 80 ) output character ( len = 80 ) string call timestamp ( ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_RULE' write ( *, '(a)' ) ' FORTRAN90 version' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Compute a Gauss-Jacobi rule for approximating' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Integral ( -1 <= x <= +1 ) (1-x)^ALPHA (1+x)^BETA f(x) dx' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' of order ORDER.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' The user specifies ORDER, ALPHA, BETA, and OUTPUT.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' OUTPUT is:' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' "C++" for printed C++ output;' write ( *, '(a)' ) ' "F77" for printed Fortran77 output;' write ( *, '(a)' ) ' "F90" for printed Fortran90 output;' write ( *, '(a)' ) ' "MAT" for printed MATLAB output;' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' or:' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' "filename" to generate 3 files:' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' filename_w.txt - the weight file' write ( *, '(a)' ) ' filename_x.txt - the abscissa file.' write ( *, '(a)' ) ' filename_r.txt - the region file.' ! ! Get the number of command line arguments. ! arg_num = iargc ( ) ! ! Get the order. ! if ( 1 <= arg_num ) then iarg = 1 call getarg ( iarg, string ) call s_to_i4 ( string, order, ierror, last ) else write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Enter the rule order ORDER:' read ( *, * ) order end if ! ! Get ALPHA ! if ( 2 <= arg_num ) then iarg = 2 call getarg ( iarg, string ) call s_to_r8 ( string, alpha, ierror, last ) else write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' ALPHA is the exponent of (1-x) in the integral:' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Integral ( -1 <= x <= +1 ) (1-x)^ALPHA (1+x)^BETA f(x) dx' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Note that -1.0 < ALPHA is required.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Enter the value of ALPHA:' read ( *, * ) alpha end if ! ! Get BETA ! if ( 3 <= arg_num ) then iarg = 3 call getarg ( iarg, string ) call s_to_r8 ( string, beta, ierror, last ) else write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' BETA is the exponent of (1+x) in the integral:' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Integral ( -1 <= x <= +1 ) (1-x)^ALPHA (1+x)^BETA f(x) dx' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Note that -1.0 < BETA is required.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Enter the value of BETA:' read ( *, * ) beta end if ! ! Get the output option or quadrature file root name: ! if ( 4 <= arg_num ) then iarg = 4 call getarg ( iarg, output ) else write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Enter OUTPUT (one of C++, F77, F90, MAT,' write ( *, '(a)' ) ' or else the "root name" of the quadrature files).' read ( *, '(a)' ) output end if ! ! Echo the input. ! write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Input summary:' write ( *, '(a)' ) ' ' write ( *, '(a,i8)' ) ' ORDER = ', order write ( *, '(a,g14.6)' ) ' ALPHA = ', alpha write ( *, '(a,g14.6)' ) ' BETA = ', beta write ( *, '(a)' ) ' OUTPUT = "' // trim ( output ) // '".' ! ! Construct the rule and output it. ! call jacobi_handle ( order, alpha, beta, output ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_RULE:' write ( *, '(a)' ) ' Normal end of execution.' write ( *, '(a)' ) ' ' call timestamp ( ) stop end subroutine ch_cap ( c ) !*****************************************************************************80 ! !! CH_CAP capitalizes a single character. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 19 July 1998 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input/output, character C, the character to capitalize. ! implicit none character c integer ( kind = 4 ) itemp itemp = ichar ( c ) if ( 97 <= itemp .and. itemp <= 122 ) then c = char ( itemp - 32 ) end if return end function ch_eqi ( c1, c2 ) !*****************************************************************************80 ! !! CH_EQI is a case insensitive comparison of two characters for equality. ! ! Example: ! ! CH_EQI ( 'A', 'a' ) is .TRUE. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 28 July 2000 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character C1, C2, the characters to compare. ! ! Output, logical CH_EQI, the result of the comparison. ! implicit none logical ch_eqi character c1 character c1_cap character c2 character c2_cap c1_cap = c1 c2_cap = c2 call ch_cap ( c1_cap ) call ch_cap ( c2_cap ) if ( c1_cap == c2_cap ) then ch_eqi = .true. else ch_eqi = .false. end if return end subroutine ch_to_digit ( c, digit ) !*****************************************************************************80 ! !! CH_TO_DIGIT returns the integer value of a base 10 digit. ! ! Example: ! ! C DIGIT ! --- ----- ! '0' 0 ! '1' 1 ! ... ... ! '9' 9 ! ' ' 0 ! 'X' -1 ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 04 August 1999 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character C, the decimal digit, '0' through '9' or blank ! are legal. ! ! Output, integer ( kind = 4 ) DIGIT, the corresponding integer value. ! If C was 'illegal', then DIGIT is -1. ! implicit none character c integer ( kind = 4 ) digit if ( lge ( c, '0' ) .and. lle ( c, '9' ) ) then digit = ichar ( c ) - 48 else if ( c == ' ' ) then digit = 0 else digit = -1 end if return end subroutine dtable_close_write ( output_unit ) !*****************************************************************************80 ! !! DTABLE_CLOSE_WRITE closes a file used to write a DTABLE. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 26 September 2006 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) OUTPUT_UNIT, the output unit that was used. ! implicit none integer ( kind = 4 ) output_unit close ( unit = output_unit ) return end subroutine dtable_data_write ( output_unit, m, n, table ) !*****************************************************************************80 ! !! DTABLE_DATA_WRITE writes DTABLE data to a file. ! ! Discussion: ! ! This routine writes a single line of output for each point, ! containing its spatial coordinates. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 07 September 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) OUTPUT_UNIT, the output unit. ! ! Input, integer ( kind = 4 ) M, the spatial dimension. ! ! Input, integer ( kind = 4 ) N, the number of points. ! ! Input, real ( kind = 8 ) TABLE(M,N), the table data. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n integer ( kind = 4 ) output_unit integer ( kind = 4 ) j character ( len = 30 ) string real ( kind = 8 ) table(m,n) ! ! Create the format string. ! write ( string, '(a1,i8,a1,i8,a1,i8,a1)' ) '(', m, 'g', 24, '.', 16, ')' call s_blank_delete ( string ) do j = 1, n write ( output_unit, string ) table(1:m,j) end do return end subroutine dtable_header_write ( output_file_name, output_unit, m, n ) !*****************************************************************************80 ! !! DTABLE_HEADER_WRITE writes the header to a DTABLE file. ! ! Discussion: ! ! The file must already be open before this routine is called. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 07 September 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) OUTPUT_FILE_NAME, the output file name. ! ! Input, integer ( kind = 4 ) OUTPUT_UNIT, the output unit. ! ! Input, integer ( kind = 4 ) M, the spatial dimension. ! ! Input, integer ( kind = 4 ) N, the number of points. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n character ( len = * ) output_file_name integer ( kind = 4 ) output_unit character ( len = 40 ) string real ( kind = 8 ), parameter :: x = 1.0D+00 call timestring ( string ) write ( output_unit, '(a)' ) '# ' // trim ( output_file_name ) write ( output_unit, '(a)' ) '# created by TABLE_IO.F90' write ( output_unit, '(a)' ) '# at ' // trim ( string ) write ( output_unit, '(a)' ) '#' write ( output_unit, '(a,i8)' ) '# Spatial dimension M = ', m write ( output_unit, '(a,i8)' ) '# Number of points N = ', n write ( output_unit, '(a,g14.6)' ) '# EPSILON (unit roundoff) = ', & epsilon ( x ) write ( output_unit, '(a)' ) '#' return end subroutine dtable_open_write ( output_file_name, output_unit ) !*****************************************************************************80 ! !! DTABLE_OPEN_WRITE opens a file to write a DTABLE. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 26 September 2006 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) OUTPUT_FILE_NAME, the output file name. ! ! Output, integer ( kind = 4 ) OUTPUT_UNIT, the output unit to be used. ! implicit none character ( len = * ) output_file_name integer ( kind = 4 ) output_status integer ( kind = 4 ) output_unit call get_unit ( output_unit ) open ( unit = output_unit, file = output_file_name, & status = 'replace', iostat = output_status ) if ( output_status /= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'DTABLE_OPEN_WRITE - Fatal error!' write ( *, '(a,i8)' ) ' Could not open the output file "' // & trim ( output_file_name ) // '" on unit ', output_unit output_unit = -1 stop end if return end subroutine dtable_write ( output_file_name, m, n, table, header ) !*****************************************************************************80 ! !! DTABLE_WRITE writes a DTABLE to a file. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 27 September 2006 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) OUTPUT_FILE_NAME, the output file name. ! ! Input, integer ( kind = 4 ) M, the spatial dimension. ! ! Input, integer ( kind = 4 ) N, the number of points. ! ! Input, real ( kind = 8 ) TABLE(M,N), the table data. ! ! Input, logical HEADER, is TRUE if the header is to be included. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n logical header character ( len = * ) output_file_name integer ( kind = 4 ) output_unit real ( kind = 8 ) table(m,n) call dtable_open_write ( output_file_name, output_unit ) if ( header ) then call dtable_header_write ( output_file_name, output_unit, m, n ) end if call dtable_data_write ( output_unit, m, n, table ) call dtable_close_write ( output_unit ) return end subroutine get_unit ( iunit ) !*****************************************************************************80 ! !! GET_UNIT returns a free FORTRAN unit number. ! ! Discussion: ! ! A "free" FORTRAN unit number is an integer between 1 and 99 which ! is not currently associated with an I/O device. A free FORTRAN unit ! number is needed in order to open a file with the OPEN command. ! ! If IUNIT = 0, then no free FORTRAN unit could be found, although ! all 99 units were checked (except for units 5, 6 and 9, which ! are commonly reserved for console I/O). ! ! Otherwise, IUNIT is an integer between 1 and 99, representing a ! free FORTRAN unit. Note that GET_UNIT assumes that units 5 and 6 ! are special, and will never return those values. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 18 September 2005 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Output, integer ( kind = 4 ) IUNIT, the free unit number. ! implicit none integer ( kind = 4 ) i integer ( kind = 4 ) ios integer ( kind = 4 ) iunit logical lopen iunit = 0 do i = 1, 99 if ( i /= 5 .and. i /= 6 .and. i /= 9 ) then inquire ( unit = i, opened = lopen, iostat = ios ) if ( ios == 0 ) then if ( .not. lopen ) then iunit = i return end if end if end if end do return end subroutine jacobi_compute ( order, alpha, beta, xtab, weight ) !*****************************************************************************80 ! !! JACOBI_COMPUTE computes a Gauss-Jacobi quadrature rule. ! ! Discussion: ! ! The integral to approximate: ! ! Integral ( -1 <= X <= 1 ) (1-X)**ALPHA * (1+X)**BETA * F(X) dX ! ! The quadrature rule: ! ! Sum ( 1 <= I <= ORDER ) WEIGHT(I) * F ( XTAB(I) ) ! ! Thanks to Xu Xiang of Fudan University for pointing out that ! an earlier implementation of this routine was incorrect! ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 14 May 2007 ! ! Author: ! ! FORTRAN77 original version by Arthur Stroud, Don Secrest ! FORTRAN90 version by John Burkardt ! ! Reference: ! ! Arthur Stroud, Don Secrest, ! Gaussian Quadrature Formulas, ! Prentice Hall, 1966, ! LC: QA299.4G3S7. ! ! Parameters: ! ! Input, integer ( kind = 4 ) ORDER, the order of the quadrature rule ! to be computed. ! ! Input, real ( kind = 8 ) ALPHA, BETA, the exponents of (1-X) and ! (1+X) in the quadrature rule. For simple Gauss-Legendre quadrature, ! set ALPHA = BETA = 0.0. -1.0 < ALPHA and -1.0 < BETA are required. ! ! Output, real ( kind = 8 ) XTAB(ORDER), the abscissas. ! ! Output, real ( kind = 8 ) WEIGHT(ORDER), the weights. ! implicit none integer ( kind = 4 ) order real ( kind = 8 ) alpha real ( kind = 8 ) an real ( kind = 8 ) b(order) real ( kind = 8 ) beta real ( kind = 8 ) bn real ( kind = 8 ) c(order) real ( kind = 8 ) cc real ( kind = 8 ) delta real ( kind = 8 ) dp2 integer ( kind = 4 ) i real ( kind = 8 ) p1 real ( kind = 8 ) r1 real ( kind = 8 ) r2 real ( kind = 8 ) r3 real ( kind = 8 ) r8_gamma real ( kind = 8 ) weight(order) real ( kind = 8 ) x real ( kind = 8 ) xtab(order) ! ! Check ALPHA and BETA. ! if ( alpha <= -1.0D+00 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_COMPUTE - Fatal error!' write ( *, '(a)' ) ' -1.0 < ALPHA is required.' stop end if if ( beta <= -1.0D+00 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_COMPUTE - Fatal error!' write ( *, '(a)' ) ' -1.0 < BETA is required.' stop end if ! ! Set the recursion coefficients. ! do i = 1, order if ( alpha + beta == 0.0D+00 .or. beta - alpha == 0.0D+00 ) then b(i) = 0.0D+00 else b(i) = ( alpha + beta ) * ( beta - alpha ) / & ( ( alpha + beta + real ( 2 * i, kind = 8 ) ) & * ( alpha + beta + real ( 2 * i - 2, kind = 8 ) ) ) end if if ( i == 1 ) then c(i) = 0.0D+00 else c(i) = 4.0D+00 * real ( i - 1, kind = 8 ) & * ( alpha + real ( i - 1, kind = 8 ) ) & * ( beta + real ( i - 1, kind = 8 ) ) & * ( alpha + beta + real ( i - 1, kind = 8 ) ) / & ( ( alpha + beta + real ( 2 * i - 1, kind = 8 ) ) & * ( alpha + beta + real ( 2 * i - 2, kind = 8 ) )**2 & * ( alpha + beta + real ( 2 * i - 3, kind = 8 ) ) ) end if end do delta = r8_gamma ( alpha + 1.0D+00 ) & * r8_gamma ( beta + 1.0D+00 ) & / r8_gamma ( alpha + beta + 2.0D+00 ) cc = delta * 2.0D+00**( alpha + beta + 1.0D+00 ) * product ( c(2:order) ) do i = 1, order if ( i == 1 ) then an = alpha / real ( order, kind = 8 ) bn = beta / real ( order, kind = 8 ) r1 = ( 1.0D+00 + alpha ) & * ( 2.78D+00 / ( 4.0D+00 + real ( order**2, kind = 8 ) ) & + 0.768D+00 * an / real ( order, kind = 8 ) ) r2 = 1.0D+00 + 1.48D+00 * an + 0.96D+00 * bn & + 0.452D+00 * an**2 + 0.83D+00 * an * bn x = ( r2 - r1 ) / r2 else if ( i == 2 ) then r1 = ( 4.1D+00 + alpha ) / & ( ( 1.0D+00 + alpha ) * ( 1.0D+00 + 0.156D+00 * alpha ) ) r2 = 1.0D+00 + 0.06D+00 * ( real ( order, kind = 8 ) - 8.0D+00 ) * & ( 1.0D+00 + 0.12D+00 * alpha ) / real ( order, kind = 8 ) r3 = 1.0D+00 + 0.012D+00 * beta * & ( 1.0D+00 + 0.25D+00 * abs ( alpha ) ) / real ( order, kind = 8 ) x = x - r1 * r2 * r3 * ( 1.0D+00 - x ) else if ( i == 3 ) then r1 = ( 1.67D+00 + 0.28D+00 * alpha ) / ( 1.0D+00 + 0.37D+00 * alpha ) r2 = 1.0D+00 + 0.22D+00 * ( real ( order, kind = 8 ) - 8.0D+00 ) & / real ( order, kind = 8 ) r3 = 1.0D+00 + 8.0D+00 * beta / & ( ( 6.28D+00 + beta ) * real ( order**2, kind = 8 ) ) x = x - r1 * r2 * r3 * ( xtab(1) - x ) else if ( i < order - 1 ) then x = 3.0D+00 * xtab(i-1) - 3.0D+00 * xtab(i-2) + xtab(i-3) else if ( i == order - 1 ) then r1 = ( 1.0D+00 + 0.235D+00 * beta ) / ( 0.766D+00 + 0.119D+00 * beta ) r2 = 1.0D+00 / ( 1.0D+00 + 0.639D+00 & * ( real ( order, kind = 8 ) - 4.0D+00 ) & / ( 1.0D+00 + 0.71D+00 * ( real ( order, kind = 8 ) - 4.0D+00 ) ) ) r3 = 1.0D+00 / ( 1.0D+00 + 20.0D+00 * alpha / ( ( 7.5D+00 + alpha ) * & real ( order**2, kind = 8 ) ) ) x = x + r1 * r2 * r3 * ( x - xtab(i-2) ) else if ( i == order ) then r1 = ( 1.0D+00 + 0.37D+00 * beta ) / ( 1.67D+00 + 0.28D+00 * beta ) r2 = 1.0D+00 / & ( 1.0D+00 + 0.22D+00 * ( real ( order, kind = 8 ) - 8.0D+00 ) & / real ( order, kind = 8 ) ) r3 = 1.0D+00 / ( 1.0D+00 + 8.0D+00 * alpha / & ( ( 6.28D+00 + alpha ) * real ( order**2, kind = 8 ) ) ) x = x + r1 * r2 * r3 * ( x - xtab(i-2) ) end if call jacobi_root ( x, order, alpha, beta, dp2, p1, b, c ) xtab(i) = x weight(i) = cc / ( dp2 * p1 ) end do ! ! Reverse the order of the data. ! xtab(1:order) = xtab(order:1:-1) weight(1:order) = weight(order:1:-1) return end subroutine jacobi_handle ( order, alpha, beta, output ) !*****************************************************************************80 ! !! JACOBI_HANDLE computes the requested rule and outputs it. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 02 September 2007 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) ORDER, the order of the rule. ! ! Input, real ( kind = 8 ) ALPHA, BETA, the parameters. ! ! Input, character ( len = * ) OUTPUT, specifies the output. ! * 'C++', print as C++ code. ! * 'F77', print as FORTRAN77 code. ! * 'F90', print as FORTRAN90 code. ! * 'MAT', print as MATLAB code. ! * file, write files 'file_w.txt', 'file_x.txt', 'file_r.txt' defining weights, ! abscissas, and region. ! implicit none real ( kind = 8 ) alpha real ( kind = 8 ) beta logical header integer ( kind = 4 ) i integer ( kind = 4 ) order character ( len = * ) output character ( len = 80 ) output_r character ( len = 80 ) output_w character ( len = 80 ) output_x real ( kind = 8 ) r(2) real ( kind = 8 ), allocatable, dimension ( : ) :: w real ( kind = 8 ), allocatable, dimension ( : ) :: x r(1) = -1.0D+00 r(2) = +1.0D+00 allocate ( w(order) ) allocate ( x(order) ) call jacobi_compute ( order, alpha, beta, x, w ) if ( output(1:3) == 'C++' .or. output(1:3) == 'c++' ) then write ( *, '(a)' ) '//' write ( *, '(a)' ) '// Weights W, abscissas X and range R' write ( *, '(a)' ) '// for a Gauss-Jacobi quadrature rule' write ( *, '(a,i8)' ) '// ORDER = ', order write ( *, '(a,g14.6)' ) '// ALPHA = ', alpha write ( *, '(a,g14.6)' ) '// BETA = ', beta write ( *, '(a)' ) '//' write ( *, '(a)' ) '// Standard rule:' write ( *, '(a)' ) '// Integral ( -1 <= x <= +1 ) (1-x)^ALPHA (1+x)^BETA f(x) dx' write ( *, '(a)' ) '// is to be approximated by' write ( *, '(a)' ) '// sum ( 1 <= I <= ORDER ) w(i) * f(x(i)).' write ( *, '(a)' ) '//' do i = 1, order write ( *, '(a,i2,a,g24.16,a)' ) ' w[', i-1, '] = ', w(i), ';' end do write ( *, '(a)' ) ' ' do i = 1, order write ( *, '(a,i2,a,g24.16,a)' ) ' x[', i-1, '] = ', x(i), ';' end do write ( *, '(a)' ) ' ' do i = 1, 2 write ( *, '(a,i2,a,g24.16,a)' ) ' r[', i-1, '] = ', r(i), ';' end do else if ( output(1:3) == 'F77' .or. output(1:3) == 'f77' ) then write ( *, '(a)' ) 'c' write ( *, '(a)' ) 'c Weights W, abscissas X and range R' write ( *, '(a)' ) 'c for a Gauss-Jacobi quadrature rule' write ( *, '(a,i8)' ) 'c ORDER = ', order write ( *, '(a,g14.6)' ) 'c ALPHA = ', alpha write ( *, '(a,g14.6)' ) 'c BETA = ', beta write ( *, '(a)' ) 'c' write ( *, '(a)' ) 'c Standard rule:' write ( *, '(a)' ) 'c Integral ( -1 <= x <= +1 ) (1-x)^ALPHA (1+x)^BETA f(x) dx' write ( *, '(a)' ) 'c is to be approximated by' write ( *, '(a)' ) 'c sum ( 1 <= I <= ORDER ) w(i) * f(x(i)).' write ( *, '(a)' ) 'c' do i = 1, order write ( *, '(a,i2,a,g24.16)' ) ' w(', i, ') = ', w(i) end do write ( *, '(a)' ) ' ' do i = 1, order write ( *, '(a,i2,a,g24.16)' ) ' x(', i, ') = ', x(i) end do write ( *, '(a)' ) ' ' do i = 1, 2 write ( *, '(a,i2,a,g24.16)' ) ' r(', i, ') = ', r(i) end do else if ( output(1:3) == 'F90' .or. output(1:3) == 'f90' ) then write ( *, '(a)' ) '!' write ( *, '(a)' ) '! Weights W, abscissas X and range R' write ( *, '(a)' ) '! for a Gauss-Jacobi quadrature rule' write ( *, '(a,i8)' ) '! ORDER = ', order write ( *, '(a,g14.6)' ) '! ALPHA = ', alpha write ( *, '(a,g14.6)' ) '! BETA = ', beta write ( *, '(a)' ) '!' write ( *, '(a)' ) '! Standard rule:' write ( *, '(a)' ) '! Integral ( -1 <= x <= +1 ) (1-x)^ALPHA (1+x)^BETA f(x) dx' write ( *, '(a)' ) '! is to be approximated by' write ( *, '(a)' ) '! sum ( 1 <= I <= ORDER ) w(i) * f(x(i)).' write ( *, '(a)' ) '!' do i = 1, order write ( *, '(a,i2,a,g24.16)' ) ' w(', i, ') = ', w(i) end do write ( *, '(a)' ) ' ' do i = 1, order write ( *, '(a,i2,a,g24.16)' ) ' x(', i, ') = ', x(i) end do write ( *, '(a)' ) ' ' do i = 1, 2 write ( *, '(a,i2,a,g24.16)' ) ' r(', i, ') = ', r(i) end do else if ( output(1:3) == 'MAT' ) then write ( *, '(a)' ) '%' write ( *, '(a)' ) '% Weights W, abscissas X and range R' write ( *, '(a)' ) '% for a Gauss-Jacobi quadrature rule' write ( *, '(a,i8)' ) '% ORDER = ', order write ( *, '(a,g14.6)' ) '% ALPHA = ', alpha write ( *, '(a,g14.6)' ) '% BETA = ', beta write ( *, '(a)' ) '%' write ( *, '(a)' ) '% Standard rule:' write ( *, '(a)' ) '% Integral ( -1 <= x <= +1 ) (1-x)^ALPHA (1+x)^BETA f(x) dx' write ( *, '(a)' ) '% is to be approximated by' write ( *, '(a)' ) '% sum ( 1 <= I <= ORDER ) w(i) * f(x(i)).' write ( *, '(a)' ) '%' do i = 1, order write ( *, '(a,i2,a,g24.16,a)' ) ' w(', i, ') = ', w(i), ';' end do write ( *, '(a)' ) ' ' do i = 1, order write ( *, '(a,i2,a,g24.16,a)' ) ' x(', i, ') = ', x(i), ';' end do write ( *, '(a)' ) ' ' do i = 1, 2 write ( *, '(a,i2,a,g24.16,a)' ) ' r(', i, ') = ', r(i), ';' end do else output_w = trim ( output ) // '_w.txt' output_x = trim ( output ) // '_x.txt' output_r = trim ( output ) // '_r.txt' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Creating quadrature files.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' "Root" file name is "' // trim ( output ) // '".' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Weight file will be "' // trim ( output_w ) // '".' write ( *, '(a)' ) ' Abscissa file will be "' // trim ( output_x ) // '".' write ( *, '(a)' ) ' Region file will be "' // trim ( output_r ) // '".' header = .false. call dtable_write ( output_w, 1, order, w, header ) call dtable_write ( output_x, 1, order, x, header ) call dtable_write ( output_r, 1, 2, r, header ) end if deallocate ( w ) deallocate ( x ) return end subroutine jacobi_recur ( p2, dp2, p1, x, order, alpha, beta, b, c ) !*****************************************************************************80 ! !! JACOBI_RECUR finds the value and derivative of a Jacobi polynomial. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 19 February 2008 ! ! Author: ! ! FORTRAN77 original version by Arthur Stroud, Don Secrest ! FORTRAN90 version by John Burkardt ! ! Reference: ! ! Arthur Stroud, Don Secrest, ! Gaussian Quadrature Formulas, ! Prentice Hall, 1966, ! LC: QA299.4G3S7. ! ! Parameters: ! ! Output, real ( kind = 8 ) P2, the value of J(ORDER)(X). ! ! Output, real ( kind = 8 ) DP2, the value of J'(ORDER)(X). ! ! Output, real ( kind = 8 ) P1, the value of J(ORDER-1)(X). ! ! Input, real ( kind = 8 ) X, the point at which polynomials are evaluated. ! ! Input, integer ( kind = 4 ) ORDER, the order of the polynomial ! to be computed. ! ! Input, real ( kind = 8 ) ALPHA, BETA, the exponents of (1+X) and ! (1-X) in the quadrature rule. ! ! Input, real ( kind = 8 ) B(ORDER), C(ORDER), the recursion ! coefficients. ! implicit none integer ( kind = 4 ) order real ( kind = 8 ) alpha real ( kind = 8 ) b(order) real ( kind = 8 ) beta real ( kind = 8 ) c(order) real ( kind = 8 ) dp0 real ( kind = 8 ) dp1 real ( kind = 8 ) dp2 integer ( kind = 4 ) i real ( kind = 8 ) p0 real ( kind = 8 ) p1 real ( kind = 8 ) p2 real ( kind = 8 ) x p1 = 1.0D+00 dp1 = 0.0D+00 p2 = x + ( alpha - beta ) / ( alpha + beta + 2.0D+00 ) dp2 = 1.0D+00 do i = 2, order p0 = p1 dp0 = dp1 p1 = p2 dp1 = dp2 p2 = ( x - b(i) ) * p1 - c(i) * p0 dp2 = ( x - b(i) ) * dp1 + p1 - c(i) * dp0 end do return end subroutine jacobi_root ( x, order, alpha, beta, dp2, p1, b, c ) !*****************************************************************************80 ! !! JACOBI_ROOT improves an approximate root of a Jacobi polynomial. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 19 February 2008 ! ! Author: ! ! FORTRAN77 original version by Arthur Stroud, Don Secrest ! FORTRAN90 version by John Burkardt ! ! Reference: ! ! Arthur Stroud, Don Secrest, ! Gaussian Quadrature Formulas, ! Prentice Hall, 1966, ! LC: QA299.4G3S7. ! ! Parameters: ! ! Input/output, real ( kind = 8 ) X, the approximate root, which ! should be improved on output. ! ! Input, integer ( kind = 4 ) ORDER, the order of the polynomial ! to be computed. ! ! Input, real ( kind = 8 ) ALPHA, BETA, the exponents of (1+X) and ! (1-X) in the quadrature rule. ! ! Output, real ( kind = 8 ) DP2, the value of J'(ORDER)(X). ! ! Output, real ( kind = 8 ) P1, the value of J(ORDER-1)(X). ! ! Input, real ( kind = 8 ) B(ORDER), C(ORDER), the recursion coefficients. ! implicit none integer ( kind = 4 ) order real ( kind = 8 ) alpha real ( kind = 8 ) b(order) real ( kind = 8 ) beta real ( kind = 8 ) c(order) real ( kind = 8 ) d real ( kind = 8 ) dp2 real ( kind = 8 ) eps real ( kind = 8 ) p1 real ( kind = 8 ) p2 real ( kind = 8 ) r8_epsilon integer ( kind = 4 ) step integer ( kind = 4 ), parameter :: step_max = 10 real ( kind = 8 ) x eps = r8_epsilon ( ) do step = 1, step_max call jacobi_recur ( p2, dp2, p1, x, order, alpha, beta, b, c ) d = p2 / dp2 x = x - d if ( abs ( d ) <= eps * ( abs ( x ) + 1.0D+00 ) ) then return end if end do return end function r8_epsilon ( ) !*****************************************************************************80 ! !! R8_EPSILON returns the R8 roundoff unit. ! ! Discussion: ! ! The roundoff unit is a number R which is a power of 2 with the ! property that, to the precision of the computer's arithmetic, ! 1 < 1 + R ! but ! 1 = ( 1 + R / 2 ) ! ! FORTRAN90 provides the superior library routine ! ! EPSILON ( X ) ! ! Modified: ! ! 15 April 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Output, real ( kind = 8 ) R8_EPSILON, the R8 round-off unit. ! implicit none real ( kind = 8 ) d real ( kind = 8 ) d_test real ( kind = 8 ) r8_epsilon d = 1.0D+00 d_test = 1.0D+00 + d / 2.0D+00 do while ( 1.0D+00 < d_test ) d = d / 2.0D+00 d_test = 1.0D+00 + d / 2.0D+00 end do r8_epsilon = d return end function r8_gamma ( x ) !*****************************************************************************80 ! !! R8_GAMMA evaluates Gamma(X) for a real argument. ! ! Discussion: ! ! This function was originally named DGAMMA. ! ! However, a number of Fortran compilers now include a library ! function of this name. To avoid conflicts, this function was ! renamed R8_GAMMA. ! ! This routine calculates the GAMMA function for a real argument X. ! Computation is based on an algorithm outlined in reference 1. ! The program uses rational functions that approximate the GAMMA ! function to at least 20 significant decimal digits. Coefficients ! for the approximation over the interval (1,2) are unpublished. ! Those for the approximation for 12 <= X are from reference 2. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! FORTRAN77 original version by William Cody, Laura Stoltz ! FORTRAN90 version by John Burkardt ! ! Reference: ! ! William Cody, ! An Overview of Software Development for Special Functions, ! in Numerical Analysis Dundee, 1975, ! edited by GA Watson, ! Lecture Notes in Mathematics 506, ! Springer, 1976. ! ! John Hart, Ward Cheney, Charles Lawson, Hans Maehly, ! Charles Mesztenyi, John Rice, Henry Thatcher, ! Christoph Witzgall, ! Computer Approximations, ! Wiley, 1968, ! LC: QA297.C64. ! ! Parameters: ! ! Input, real ( kind = 8 ) X, the argument of the function. ! ! Output, real ( kind = 8 ) R8_GAMMA, the value of the function. ! implicit none ! ! Coefficients for minimax approximation over (12, INF). ! real ( kind = 8 ), dimension ( 7 ) :: c = (/ & -1.910444077728D-03, & 8.4171387781295D-04, & -5.952379913043012D-04, & 7.93650793500350248D-04, & -2.777777777777681622553D-03, & 8.333333333333333331554247D-02, & 5.7083835261D-03 /) real ( kind = 8 ), parameter :: eps = 2.22D-16 real ( kind = 8 ) fact real ( kind = 8 ), parameter :: half = 0.5D+00 integer ( kind = 4 ) i integer ( kind = 4 ) n real ( kind = 8 ), parameter :: one = 1.0D+00 real ( kind = 8 ), dimension ( 8 ) :: p = (/ & -1.71618513886549492533811D+00, & 2.47656508055759199108314D+01, & -3.79804256470945635097577D+02, & 6.29331155312818442661052D+02, & 8.66966202790413211295064D+02, & -3.14512729688483675254357D+04, & -3.61444134186911729807069D+04, & 6.64561438202405440627855D+04 /) logical parity real ( kind = 8 ), parameter :: pi = 3.1415926535897932384626434D+00 real ( kind = 8 ), dimension ( 8 ) :: q = (/ & -3.08402300119738975254353D+01, & 3.15350626979604161529144D+02, & -1.01515636749021914166146D+03, & -3.10777167157231109440444D+03, & 2.25381184209801510330112D+04, & 4.75584627752788110767815D+03, & -1.34659959864969306392456D+05, & -1.15132259675553483497211D+05 /) real ( kind = 8 ) r8_gamma real ( kind = 8 ) res real ( kind = 8 ), parameter :: sqrtpi = 0.9189385332046727417803297D+00 real ( kind = 8 ) sum real ( kind = 8 ), parameter :: twelve = 12.0D+00 real ( kind = 8 ), parameter :: two = 2.0D+00 real ( kind = 8 ) x real ( kind = 8 ), parameter :: xbig = 171.624D+00 real ( kind = 8 ) xden real ( kind = 8 ), parameter :: xinf = 1.0D+30 real ( kind = 8 ), parameter :: xminin = 2.23D-308 real ( kind = 8 ) xnum real ( kind = 8 ) y real ( kind = 8 ) y1 real ( kind = 8 ) ysq real ( kind = 8 ) z real ( kind = 8 ), parameter :: zero = 0.0D+00 parity = .false. fact = one n = 0 y = x ! ! Argument is negative. ! if ( y <= zero ) then y = - x y1 = aint ( y ) res = y - y1 if ( res /= zero ) then if ( y1 /= aint ( y1 * half ) * two ) then parity = .true. end if fact = - pi / sin ( pi * res ) y = y + one else res = xinf r8_gamma = res return end if end if ! ! Argument is positive. ! if ( y < eps ) then ! ! Argument < EPS. ! if ( xminin <= y ) then res = one / y else res = xinf r8_gamma = res return end if else if ( y < twelve ) then y1 = y ! ! 0.0 < argument < 1.0. ! if ( y < one ) then z = y y = y + one ! ! 1.0 < argument < 12.0. ! Reduce argument if necessary. ! else n = int ( y ) - 1 y = y - real ( n, kind = 8 ) z = y - one end if ! ! Evaluate approximation for 1.0 < argument < 2.0. ! xnum = zero xden = one do i = 1, 8 xnum = ( xnum + p(i) ) * z xden = xden * z + q(i) end do res = xnum / xden + one ! ! Adjust result for case 0.0 < argument < 1.0. ! if ( y1 < y ) then res = res / y1 ! ! Adjust result for case 2.0 < argument < 12.0. ! else if ( y < y1 ) then do i = 1, n res = res * y y = y + one end do end if else ! ! Evaluate for 12.0 <= argument. ! if ( y <= xbig ) then ysq = y * y sum = c(7) do i = 1, 6 sum = sum / ysq + c(i) end do sum = sum / y - y + sqrtpi sum = sum + ( y - half ) * log ( y ) res = exp ( sum ) else res = xinf r8_gamma = res return end if end if ! ! Final adjustments and return. ! if ( parity ) then res = - res end if if ( fact /= one ) then res = fact / res end if r8_gamma = res return end subroutine s_blank_delete ( s ) !*****************************************************************************80 ! !! S_BLANK_DELETE removes blanks from a string, left justifying the remainder. ! ! Discussion: ! ! All TAB characters are also removed. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 26 July 1998 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input/output, character ( len = * ) S, the string to be transformed. ! implicit none character c integer ( kind = 4 ) get integer ( kind = 4 ) put integer ( kind = 4 ) nchar character ( len = * ) s character, parameter :: TAB = char ( 9 ) put = 0 nchar = len_trim ( s ) do get = 1, nchar c = s(get:get) if ( c /= ' ' .and. c /= TAB ) then put = put + 1 s(put:put) = c end if end do s(put+1:nchar) = ' ' return end subroutine s_to_i4 ( s, ival, ierror, length ) !*****************************************************************************80 ! !! S_TO_I4 reads an I4 from a string. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 28 June 2000 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) S, a string to be examined. ! ! Output, integer ( kind = 4 ) IVAL, the integer value read from the string. ! If the string is blank, then IVAL will be returned 0. ! ! Output, integer ( kind = 4 ) IERROR, an error flag. ! 0, no error. ! 1, an error occurred. ! ! Output, integer ( kind = 4 ) LENGTH, the number of characters of S ! used to make IVAL. ! implicit none character c integer ( kind = 4 ) i integer ( kind = 4 ) ierror integer ( kind = 4 ) isgn integer ( kind = 4 ) istate integer ( kind = 4 ) ival integer ( kind = 4 ) length character ( len = * ) s ierror = 0 istate = 0 isgn = 1 ival = 0 do i = 1, len_trim ( s ) c = s(i:i) ! ! Haven't read anything. ! if ( istate == 0 ) then if ( c == ' ' ) then else if ( c == '-' ) then istate = 1 isgn = -1 else if ( c == '+' ) then istate = 1 isgn = + 1 else if ( lle ( '0', c ) .and. lle ( c, '9' ) ) then istate = 2 ival = ichar ( c ) - ichar ( '0' ) else ierror = 1 return end if ! ! Have read the sign, expecting digits. ! else if ( istate == 1 ) then if ( c == ' ' ) then else if ( lle ( '0', c ) .and. lle ( c, '9' ) ) then istate = 2 ival = ichar ( c ) - ichar ( '0' ) else ierror = 1 return end if ! ! Have read at least one digit, expecting more. ! else if ( istate == 2 ) then if ( lle ( '0', c ) .and. lle ( c, '9' ) ) then ival = 10 * ival + ichar ( c ) - ichar ( '0' ) else ival = isgn * ival length = i - 1 return end if end if end do ! ! If we read all the characters in the string, see if we're OK. ! if ( istate == 2 ) then ival = isgn * ival length = len_trim ( s ) else ierror = 1 length = 0 end if return end subroutine s_to_r8 ( s, dval, ierror, length ) !*****************************************************************************80 ! !! S_TO_R8 reads an R8 from a string. ! ! Discussion: ! ! The routine will read as many characters as possible until it reaches ! the end of the string, or encounters a character which cannot be ! part of the number. ! ! Legal input is: ! ! 1 blanks, ! 2 '+' or '-' sign, ! 2.5 blanks ! 3 integer part, ! 4 decimal point, ! 5 fraction part, ! 6 'E' or 'e' or 'D' or 'd', exponent marker, ! 7 exponent sign, ! 8 exponent integer part, ! 9 exponent decimal point, ! 10 exponent fraction part, ! 11 blanks, ! 12 final comma or semicolon, ! ! with most quantities optional. ! ! Example: ! ! S DVAL ! ! '1' 1.0 ! ' 1 ' 1.0 ! '1A' 1.0 ! '12,34,56' 12.0 ! ' 34 7' 34.0 ! '-1E2ABCD' -100.0 ! '-1X2ABCD' -1.0 ! ' 2E-1' 0.2 ! '23.45' 23.45 ! '-4.2E+2' -420.0 ! '17d2' 1700.0 ! '-14e-2' -0.14 ! 'e2' 100.0 ! '-12.73e-9.23' -12.73 * 10.0**(-9.23) ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 07 September 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) S, the string containing the ! data to be read. Reading will begin at position 1 and ! terminate at the end of the string, or when no more ! characters can be read to form a legal real. Blanks, ! commas, or other nonnumeric data will, in particular, ! cause the conversion to halt. ! ! Output, real ( kind = 8 ) DVAL, the value read from the string. ! ! Output, integer ( kind = 4 ) IERROR, error flag. ! 0, no errors occurred. ! 1, 2, 6 or 7, the input number was garbled. The ! value of IERROR is the last type of input successfully ! read. For instance, 1 means initial blanks, 2 means ! a plus or minus sign, and so on. ! ! Output, integer ( kind = 4 ) LENGTH, the number of characters read ! to form the number, including any terminating ! characters such as a trailing comma or blanks. ! implicit none logical ch_eqi character c real ( kind = 8 ) dval integer ( kind = 4 ) ierror integer ( kind = 4 ) ihave integer ( kind = 4 ) isgn integer ( kind = 4 ) iterm integer ( kind = 4 ) jbot integer ( kind = 4 ) jsgn integer ( kind = 4 ) jtop integer ( kind = 4 ) length integer ( kind = 4 ) nchar integer ( kind = 4 ) ndig real ( kind = 8 ) rbot real ( kind = 8 ) rexp real ( kind = 8 ) rtop character ( len = * ) s nchar = len_trim ( s ) ierror = 0 dval = 0.0D+00 length = -1 isgn = 1 rtop = 0 rbot = 1 jsgn = 1 jtop = 0 jbot = 1 ihave = 1 iterm = 0 do length = length + 1 if ( nchar < length+1 ) then exit end if c = s(length+1:length+1) ! ! Blank character. ! if ( c == ' ' ) then if ( ihave == 2 ) then else if ( ihave == 6 .or. ihave == 7 ) then iterm = 1 else if ( 1 < ihave ) then ihave = 11 end if ! ! Comma. ! else if ( c == ',' .or. c == ';' ) then if ( ihave /= 1 ) then iterm = 1 ihave = 12 length = length + 1 end if ! ! Minus sign. ! else if ( c == '-' ) then if ( ihave == 1 ) then ihave = 2 isgn = -1 else if ( ihave == 6 ) then ihave = 7 jsgn = -1 else iterm = 1 end if ! ! Plus sign. ! else if ( c == '+' ) then if ( ihave == 1 ) then ihave = 2 else if ( ihave == 6 ) then ihave = 7 else iterm = 1 end if ! ! Decimal point. ! else if ( c == '.' ) then if ( ihave < 4 ) then ihave = 4 else if ( 6 <= ihave .and. ihave <= 8 ) then ihave = 9 else iterm = 1 end if ! ! Scientific notation exponent marker. ! else if ( ch_eqi ( c, 'E' ) .or. ch_eqi ( c, 'D' ) ) then if ( ihave < 6 ) then ihave = 6 else iterm = 1 end if ! ! Digit. ! else if ( ihave < 11 .and. lle ( '0', c ) .and. lle ( c, '9' ) ) then if ( ihave <= 2 ) then ihave = 3 else if ( ihave == 4 ) then ihave = 5 else if ( ihave == 6 .or. ihave == 7 ) then ihave = 8 else if ( ihave == 9 ) then ihave = 10 end if call ch_to_digit ( c, ndig ) if ( ihave == 3 ) then rtop = 10.0D+00 * rtop + real ( ndig, kind = 8 ) else if ( ihave == 5 ) then rtop = 10.0D+00 * rtop + real ( ndig, kind = 8 ) rbot = 10.0D+00 * rbot else if ( ihave == 8 ) then jtop = 10 * jtop + ndig else if ( ihave == 10 ) then jtop = 10 * jtop + ndig jbot = 10 * jbot end if ! ! Anything else is regarded as a terminator. ! else iterm = 1 end if ! ! If we haven't seen a terminator, and we haven't examined the ! entire string, go get the next character. ! if ( iterm == 1 ) then exit end if end do ! ! If we haven't seen a terminator, and we have examined the ! entire string, then we're done, and LENGTH is equal to NCHAR. ! if ( iterm /= 1 .and. length+1 == nchar ) then length = nchar end if ! ! Number seems to have terminated. Have we got a legal number? ! Not if we terminated in states 1, 2, 6 or 7! ! if ( ihave == 1 .or. ihave == 2 .or. ihave == 6 .or. ihave == 7 ) then ierror = ihave write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'S_TO_R8 - Serious error!' write ( *, '(a)' ) ' Illegal or nonnumeric input:' write ( *, '(a)' ) ' ' // trim ( s ) return end if ! ! Number seems OK. Form it. ! if ( jtop == 0 ) then rexp = 1.0D+00 else if ( jbot == 1 ) then rexp = 10.0D+00 ** ( jsgn * jtop ) else rexp = 10.0D+00 ** ( real ( jsgn * jtop, kind = 8 ) & / real ( jbot, kind = 8 ) ) end if end if dval = real ( isgn, kind = 8 ) * rexp * rtop / rbot return end subroutine timestamp ( ) !*****************************************************************************80 ! !! TIMESTAMP prints the current YMDHMS date as a time stamp. ! ! Example: ! ! May 31 2001 9:45:54.872 AM ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 31 May 2001 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! None ! implicit none character ( len = 8 ) ampm integer ( kind = 4 ) d character ( len = 8 ) date integer ( kind = 4 ) h integer ( kind = 4 ) m integer ( kind = 4 ) mm character ( len = 9 ), parameter, dimension(12) :: month = (/ & 'January ', 'February ', 'March ', 'April ', & 'May ', 'June ', 'July ', 'August ', & 'September', 'October ', 'November ', 'December ' /) integer ( kind = 4 ) n integer ( kind = 4 ) s character ( len = 10 ) time integer ( kind = 4 ) values(8) integer ( kind = 4 ) y character ( len = 5 ) zone call date_and_time ( date, time, zone, 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 ( *, '(a,1x,i2,1x,i4,2x,i2,a1,i2.2,a1,i2.2,a1,i3.3,1x,a)' ) & trim ( month(m) ), d, y, h, ':', n, ':', s, '.', mm, trim ( ampm ) return end subroutine timestring ( string ) !*****************************************************************************80 ! !! TIMESTRING writes the current YMDHMS date into a string. ! ! Example: ! ! STRING = 'May 31 2001 9:45:54.872 AM' ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 15 March 2003 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Output, character ( len = * ) STRING, contains the date information. ! A character length of 40 should always be sufficient. ! implicit none character ( len = 8 ) ampm integer ( kind = 4 ) d character ( len = 8 ) date integer ( kind = 4 ) h integer ( kind = 4 ) m integer ( kind = 4 ) mm character ( len = 9 ), parameter, dimension(12) :: month = (/ & 'January ', 'February ', 'March ', 'April ', & 'May ', 'June ', 'July ', 'August ', & 'September', 'October ', 'November ', 'December ' /) integer ( kind = 4 ) n integer ( kind = 4 ) s character ( len = * ) string character ( len = 10 ) time integer ( kind = 4 ) values(8) integer ( kind = 4 ) y character ( len = 5 ) zone call date_and_time ( date, time, zone, 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 ( string, '(a,1x,i2,1x,i4,2x,i2,a1,i2.2,a1,i2.2,a1,i3.3,1x,a)' ) & trim ( month(m) ), d, y, h, ':', n, ':', s, '.', mm, trim ( ampm ) return end