program main !*****************************************************************************80 ! !! jacobi_exactness() investigates a standard Gauss-Jacobi quadrature rule. ! ! Discussion: ! ! This program investigates a standard Gauss-Jacobi quadrature rule ! by using it to integrate monomials over [-1,+1], and comparing the ! approximate result to the known exact value. ! ! The user specifies: ! * the "root" name of the R, W and X files that specify the rule; ! * DEGREE_MAX, the maximum monomial degree to be checked. ! * ALPHA, the exponent of the (1-X) factor; ! * BETA, the exponent of the (1+X) factor. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 05 August 2009 ! ! Author: ! ! John Burkardt ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) real ( kind = rk ) alpha integer arg_num real ( kind = rk ) beta integer degree integer degree_max integer dim_num integer dim_num2 integer expon integer i integer iarg integer iargc integer ierror integer last integer order integer point_num real ( kind = rk ) quad_error character ( len = 255 ) quad_filename character ( len = 255 ) quad_r_filename character ( len = 255 ) quad_w_filename character ( len = 255 ) quad_x_filename real ( kind = rk ), allocatable, dimension ( : ) :: r character ( len = 255 ) string real ( kind = rk ), allocatable, dimension ( : ) :: w real ( kind = rk ), allocatable, dimension ( : ) :: x call timestamp ( ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'jacobi_exactness():' write ( *, '(a)' ) ' FORTRAN90 version' write ( *, '(a)' ) ' Investigate the polynomial exactness of a Gauss-Jacobi' write ( *, '(a)' ) ' quadrature rule by integrating weighted ' write ( *, '(a)' ) & ' monomials up to a given degree over the [-1,+1] interval.' ! ! Get the number of command line arguments. ! arg_num = iargc ( ) ! ! Get the quadrature file root name: ! if ( 1 <= arg_num ) then iarg = 1 call getarg ( iarg, quad_filename ) else write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_EXACTNESS:' write ( *, '(a)' ) ' Enter the "root" name of the quadrature files.' read ( *, '(a)' ) quad_filename end if ! ! Create the names of: ! the quadrature X file; ! the quadrature W file; ! the quadrature R file; ! quad_x_filename = trim ( quad_filename ) // '_x.txt' quad_w_filename = trim ( quad_filename ) // '_w.txt' quad_r_filename = trim ( quad_filename ) // '_r.txt' ! ! The second command line argument is the maximum degree. ! if ( 2 <= arg_num ) then iarg = 2 call getarg ( iarg, string ) call s_to_i4 ( string, degree_max, ierror, last ) else write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_EXACTNESS:' write ( *, '(a)' ) ' Please enter the maximum degree to check.' read ( *, * ) degree_max end if ! ! The third command line argument is ALPHA. ! if ( 3 <= arg_num ) then iarg = 3 call getarg ( iarg, string ) call s_to_r8 ( string, alpha, ierror, last ) else write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_EXACTNESS:' write ( *, '(a)' ) & ' ALPHA is the power of (1-X) in the weighting function.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' ALPHA is a real number greater than -1.0.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Please enter ALPHA:' read ( *, * ) alpha end if ! ! The fourth command line argument is BETA. ! if ( 4 <= arg_num ) then iarg = 4 call getarg ( iarg, string ) call s_to_r8 ( string, beta, ierror, last ) else write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_EXACTNESS:' write ( *, '(a)' ) ' BETA is the power of (1+X) in the weighting function.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' BETA is a real number greater than -1.0.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Please enter BETA:' read ( *, * ) beta end if ! ! Summarize the input. ! write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_EXACTNESS: User input:' write ( *, '(a)' ) ' Quadrature rule X file = "' & // trim ( quad_x_filename ) // '".' write ( *, '(a)' ) ' Quadrature rule W file = "' & // trim ( quad_w_filename ) // '".' write ( *, '(a)' ) ' Quadrature rule R file = "' & // trim ( quad_r_filename ) // '".' write ( *, '(a,i8)' ) ' Maximum degree to check = ', & degree_max write ( *, '(a,g14.6)' ) ' Exponent of (1-x), ALPHA = ', alpha write ( *, '(a,g14.6)' ) ' Exponent of (1+x), BETA = ', beta ! ! Read the X file. ! call r8mat_header_read ( quad_x_filename, dim_num, order ) if ( dim_num /= 1 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_EXACTNESS - Fatal error!' write ( *, '(a)' ) ' The spatial dimension should be 1.' write ( *, '(a,i8)' ) & ' The implicit input dimension was DIM_NUM = ', dim_num stop 1 end if write ( *, '(a)' ) ' ' write ( *, '(a,i8)' ) ' Spatial dimension = ', dim_num write ( *, '(a,i8)' ) ' Number of points = ', order allocate ( x(order) ) call r8mat_data_read ( quad_x_filename, dim_num, order, x ) ! ! Read the W file. ! call r8mat_header_read ( quad_w_filename, dim_num2, point_num ) if ( dim_num2 /= 1 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_EXACTNESS - Fatal error!' write ( *, '(a)' ) ' The quadrature weight file should have exactly' write ( *, '(a)' ) ' one value on each line.' stop 1 end if if ( point_num /= order ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_EXACTNESS - Fatal error!' write ( *, '(a)' ) ' The quadrature weight file should have exactly' write ( *, '(a)' ) ' the same number of lines as the abscissa file.' stop 1 end if allocate ( w(order) ) call r8mat_data_read ( quad_w_filename, 1, order, w ) ! ! Read the R file. ! call r8mat_header_read ( quad_r_filename, dim_num2, point_num ) if ( dim_num2 /= dim_num ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_EXACTNESS - Fatal error!' write ( *, '(a)' ) ' The quadrature region file should have the' write ( *, '(a)' ) ' same number of values on each line as the' write ( *, '(a)' ) ' abscissa file does.' stop 1 end if if ( point_num /= 2 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_EXACTNESS - Fatal error!' write ( *, '(a)' ) ' The quadrature region file should have two lines.' stop 1 end if allocate ( r(2) ) call r8mat_data_read ( quad_r_filename, dim_num, 2, r ) ! ! Print the input quadrature rule. ! write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' The quadrature rule to be tested is' write ( *, '(a)' ) ' a Gauss-Jacobi 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)' ) ' ' write ( *, '(a)' ) ' Weights W:' write ( *, '(a)' ) ' ' do i = 1, order write ( *, '(a,i2,a,g24.16)' ) ' w(', i, ') = ', w(i) end do write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Abscissas X:' write ( *, '(a)' ) ' ' do i = 1, order write ( *, '(a,i2,a,g24.16)' ) ' x(', i, ') = ', x(i) end do write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Region R:' write ( *, '(a)' ) ' ' do i = 1, 2 write ( *, '(a,i2,a,g24.16)' ) ' r(', i, ') = ', r(i) end do ! ! Explore the monomials. ! write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' A Gauss-Jacobi rule would be able to exactly' write ( *, '(a,i8)' ) & ' integrate monomials up to and including degree = ', 2 * order - 1 write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Error Degree Exponents' write ( *, '(a)' ) ' ' do degree = 0, degree_max expon = degree call monomial_quadrature_jacobi ( expon, alpha, beta, order, w, x, & quad_error ) write ( *, '(2x,f24.16,3x,i2,4x,i2)' ) & quad_error, degree, expon end do ! ! Free memory. ! deallocate ( w ) deallocate ( x ) deallocate ( r ) ! ! Terminate. ! write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'JACOBI_EXACTNESS:' write ( *, '(a)' ) ' Normal end of execution.' write ( *, '(a)' ) ' ' call timestamp ( ) stop 0 end subroutine ch_cap ( c ) !*****************************************************************************80 ! !! CH_CAP capitalizes a single character. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input/output, character C, the character to capitalize. ! implicit none character c integer 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 MIT license. ! ! Modified: ! ! 11 February 2008 ! ! 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 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 MIT license. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character C, the decimal digit, '0' through '9' or blank ! are legal. ! ! Output, integer DIGIT, the corresponding value. ! If C was 'illegal', then DIGIT is -1. ! implicit none character c integer 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 file_column_count ( input_file_name, column_num ) !*****************************************************************************80 ! !! FILE_COLUMN_COUNT counts the number of columns in the first line of a file. ! ! Discussion: ! ! The file is assumed to be a simple text file. ! ! Most lines of the file is presumed to consist of COLUMN_NUM words, ! separated by spaces. There may also be some blank lines, and some ! comment lines, ! which have a "#" in column 1. ! ! The routine tries to find the first non-comment non-blank line and ! counts the number of words in that line. ! ! If all lines are blanks or comments, it goes back and tries to analyze ! a comment line. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILE_NAME, the name of the file. ! ! Output, integer COLUMN_NUM, the number of columns in the file. ! implicit none integer column_num logical got_one character ( len = * ) input_file_name integer input_status integer input_unit character ( len = 255 ) line ! ! Open the file. ! call get_unit ( input_unit ) open ( unit = input_unit, file = input_file_name, status = 'old', & form = 'formatted', access = 'sequential', iostat = input_status ) if ( input_status /= 0 ) then column_num = -1 write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'FILE_COLUMN_COUNT - Fatal error!' write ( *, '(a,i8)' ) ' Could not open the input file "' & // trim ( input_file_name ) // '" on unit ', input_unit return end if ! ! Read one line, but skip blank lines and comment lines. ! got_one = .false. do read ( input_unit, '(a)', iostat = input_status ) line if ( input_status /= 0 ) then exit end if if ( len_trim ( line ) == 0 ) then cycle end if if ( line(1:1) == '#' ) then cycle end if got_one = .true. exit end do if ( .not. got_one ) then rewind ( input_unit ) do read ( input_unit, '(a)', iostat = input_status ) line if ( input_status /= 0 ) then exit end if if ( len_trim ( line ) == 0 ) then cycle end if got_one = .true. exit end do end if close ( unit = input_unit ) if ( .not. got_one ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'FILE_COLUMN_COUNT - Warning!' write ( *, '(a)' ) ' The file does not seem to contain any data.' column_num = -1 return end if call s_word_count ( line, column_num ) return end subroutine file_row_count ( input_file_name, row_num ) !*****************************************************************************80 ! !! FILE_ROW_COUNT counts the number of row records in a file. ! ! Discussion: ! ! It does not count lines that are blank, or that begin with a ! comment symbol '#'. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILE_NAME, the name of the input file. ! ! Output, integer ROW_NUM, the number of rows found. ! implicit none integer bad_num integer comment_num integer ierror character ( len = * ) input_file_name integer input_status integer input_unit character ( len = 255 ) line integer record_num integer row_num call get_unit ( input_unit ) open ( unit = input_unit, file = input_file_name, status = 'old', & iostat = input_status ) if ( input_status /= 0 ) then row_num = -1; ierror = 1 write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'FILE_ROW_COUNT - Fatal error!' write ( *, '(a,i8)' ) ' Could not open the input file "' // & trim ( input_file_name ) // '" on unit ', input_unit stop end if comment_num = 0 row_num = 0 record_num = 0 bad_num = 0 do read ( input_unit, '(a)', iostat = input_status ) line if ( input_status /= 0 ) then ierror = record_num exit end if record_num = record_num + 1 if ( line(1:1) == '#' ) then comment_num = comment_num + 1 cycle end if if ( len_trim ( line ) == 0 ) then comment_num = comment_num + 1 cycle end if row_num = row_num + 1 end do close ( unit = input_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 MIT license. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Output, integer IUNIT, the free unit number. ! implicit none integer i integer ios integer 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_integral ( expon, alpha, beta, value ) !*****************************************************************************80 ! !! JACOBI_INTEGRAL evaluates the integral of a monomial with Jacobi weight. ! ! Discussion: ! ! VALUE = Integral ( -1 <= X <= +1 ) x^EXPON (1-x)^ALPHA (1+x)^BETA dx ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer EXPON, the exponent. ! ! Input, real ( kind = rk ) ALPHA, the exponent of (1-X) in the weight factor. ! ! Input, real ( kind = rk ) BETA, the exponent of (1+X) in the weight factor. ! ! Output, real ( kind = rk ) VALUE, the value of the integral. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) real ( kind = rk ) alpha real ( kind = rk ) arg1 real ( kind = rk ) arg2 real ( kind = rk ) arg3 real ( kind = rk ) arg4 real ( kind = rk ) beta real ( kind = rk ) c integer expon real ( kind = rk ) s real ( kind = rk ) value real ( kind = rk ) value1 real ( kind = rk ) value2 c = real ( expon, kind = rk ) if ( mod ( expon, 2 ) == 0 ) then s = +1.0D+00 else s = -1.0D+00 end if arg1 = - alpha arg2 = 1.0D+00 + c arg3 = 2.0D+00 + beta + c arg4 = - 1.0D+00 call r8_hyper_2f1 ( arg1, arg2, arg3, arg4, value1 ) arg1 = - beta arg2 = 1.0D+00 + c arg3 = 2.0D+00 + alpha + c arg4 = - 1.0D+00 call r8_hyper_2f1 ( arg1, arg2, arg3, arg4, value2 ) value = gamma ( 1.0D+00 + c ) * ( & s * gamma ( 1.0D+00 + beta ) * value1 & / gamma ( 2.0D+00 + beta + c ) & + gamma ( 1.0D+00 + alpha ) * value2 & / gamma ( 2.0D+00 + alpha + c ) ) return end subroutine monomial_quadrature_jacobi ( expon, alpha, beta, order, w, x, & quad_error ) !*****************************************************************************80 ! !! MONOMIAL_QUADRATURE_JACOBI applies a quadrature rule to a monomial. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer EXPON, the exponent. ! ! Input, real ( kind = rk ) ALPHA, the exponent of (1-X) in the weight factor. ! ! Input, real ( kind = rk ) BETA, the exponent of (1+X) in the weight factor. ! ! Input, integer ORDER, the number of points in the rule. ! ! Input, real ( kind = rk ) W(ORDER), the quadrature weights. ! ! Input, real ( kind = rk ) X(ORDER), the quadrature points. ! ! Output, real ( kind = rk ) QUAD_ERROR, the quadrature error. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) integer order real ( kind = rk ) alpha real ( kind = rk ) beta real ( kind = rk ) exact integer expon real ( kind = rk ) quad real ( kind = rk ) quad_error real ( kind = rk ) value(order) real ( kind = rk ) w(order) real ( kind = rk ) x(order) ! ! Get the exact value of the integral. ! call jacobi_integral ( expon, alpha, beta, exact ) ! ! Evaluate the monomial at the quadrature points. ! value(1:order) = x(1:order)**expon ! ! Compute the weighted sum. ! quad = dot_product ( w, value ) ! ! Error: ! if ( exact == 0.0D+00 ) then quad_error = abs ( quad - exact ) else quad_error = abs ( ( quad - exact ) / exact ) end if return end subroutine r8_hyper_2f1 ( a_input, b_input, c_input, x_input, hf ) !*****************************************************************************80 ! !! R8_HYPER_2F1 evaluates the hypergeometric function F(A,B,C,X). ! ! Discussion: ! ! A minor bug was corrected. The HW variable, used in several places as ! the "old" value of a quantity being iteratively improved, was not ! being initialized. JVB, 11 February 2008. ! ! The original version of this program allowed the input arguments to ! be modified, although they were restored to their input values before exit. ! This is unacceptable if the input arguments are allowed to be constants. ! The code has been modified so that the input arguments are never modified. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 11 October 2008 ! ! Author: ! ! Original FORTRAN77 version by Shanjie Zhang, Jianming Jin. ! FORTRAN90 version by John Burkardt. ! ! The original FORTRAN77 version of this routine is copyrighted by ! Shanjie Zhang and Jianming Jin. However, they give permission to ! incorporate this routine into a user program provided that the copyright ! is acknowledged. ! ! Reference: ! ! Shanjie Zhang, Jianming Jin, ! Computation of Special Functions, ! Wiley, 1996, ! ISBN: 0-471-11963-6, ! LC: QA351.C45 ! ! Parameters: ! ! Input, real ( kind = rk ) A_INPUT, B_INPUT, C_INPUT, X_INPUT, ! the arguments of the function. The user is allowed to pass these ! values as constants or variables. ! C_INPUT must not be equal to a nonpositive integer. ! X_INPUT < 1. ! ! Output, real ( kind = rk ) HF, the value of the function. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) real ( kind = rk ) a real ( kind = rk ) a_input real ( kind = rk ) a0 real ( kind = rk ) aa real ( kind = rk ) b real ( kind = rk ) b_input real ( kind = rk ) bb real ( kind = rk ) c real ( kind = rk ) c_input real ( kind = rk ) c0 real ( kind = rk ) c1 real ( kind = rk ), parameter :: el = 0.5772156649015329D+00 real ( kind = rk ) eps real ( kind = rk ) f0 real ( kind = rk ) f1 real ( kind = rk ) g0 real ( kind = rk ) g1 real ( kind = rk ) g2 real ( kind = rk ) g3 real ( kind = rk ) ga real ( kind = rk ) gabc real ( kind = rk ) gam real ( kind = rk ) gb real ( kind = rk ) gbm real ( kind = rk ) gc real ( kind = rk ) gca real ( kind = rk ) gcab real ( kind = rk ) gcb real ( kind = rk ) gm real ( kind = rk ) hf real ( kind = rk ) hw integer j integer k logical l0 logical l1 logical l2 logical l3 logical l4 logical l5 integer m integer nm real ( kind = rk ) pa real ( kind = rk ) pb real ( kind = rk ), parameter :: pi = 3.141592653589793D+00 real ( kind = rk ) r real ( kind = rk ) r0 real ( kind = rk ) r1 real ( kind = rk ) r8_psi real ( kind = rk ) rm real ( kind = rk ) rp real ( kind = rk ) sm real ( kind = rk ) sp real ( kind = rk ) sp0 real ( kind = rk ) x real ( kind = rk ) x_input real ( kind = rk ) x1 ! ! Immediately copy the input arguments! ! a = a_input b = b_input c = c_input x = x_input l0 = ( c == aint ( c ) ) .and. ( c < 0.0D+00 ) l1 = ( 1.0D+00 - x < 1.0D-15 ) .and. ( c - a - b <= 0.0D+00 ) l2 = ( a == aint ( a ) ) .and. ( a < 0.0D+00 ) l3 = ( b == aint ( b ) ) .and. ( b < 0.0D+00 ) l4 = ( c - a == aint ( c - a ) ) .and. ( c - a <= 0.0D+00 ) l5 = ( c - b == aint ( c - b ) ) .and. ( c - b <= 0.0D+00 ) if ( l0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8_HYPER_2F1 - Fatal error!' write ( *, '(a)' ) ' Integral C < 0.' write ( *, '(a)' ) ' The hypergeometric series is divergent.' stop end if if ( l1 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8_HYPER_2F1 - Fatal error!' write ( *, '(a)' ) ' The hypergeometric series is divergent.' write ( *, '(a)' ) ' 1 - X < 0, C - A - B < 0.' stop end if if ( 0.95D+00 < x ) then eps = 1.0D-08 else eps = 1.0D-15 end if if ( x == 0.0D+00 .or. a == 0.0D+00 .or. b == 0.0D+00 ) then hf = 1.0D+00 return else if ( 1.0D+00 - x == eps .and. 0.0D+00 < c - a - b ) then gc = gamma ( c ) gcab = gamma ( c - a - b ) gca = gamma ( c - a ) gcb = gamma ( c - b ) hf = gc * gcab / ( gca * gcb ) return else if ( 1.0D+00 + x <= eps .and. abs ( c - a + b - 1.0D+00 ) <= eps ) then g0 = sqrt ( pi ) * 2.0D+00**( - a ) g1 = gamma ( c ) g2 = gamma ( 1.0D+00 + a / 2.0D+00 - b ) g3 = gamma ( 0.5D+00 + 0.5D+00 * a ) hf = g0 * g1 / ( g2 * g3 ) return else if ( l2 .or. l3 ) then if ( l2 ) then nm = int ( abs ( a ) ) end if if ( l3 ) then nm = int ( abs ( b ) ) end if hf = 1.0D+00 r = 1.0D+00 do k = 1, nm r = r * ( a + k - 1.0D+00 ) * ( b + k - 1.0D+00 ) & / ( k * ( c + k - 1.0D+00 ) ) * x hf = hf + r end do return else if ( l4 .or. l5 ) then if ( l4 ) then nm = int ( abs ( c - a ) ) end if if ( l5 ) then nm = int ( abs ( c - b ) ) end if hf = 1.0D+00 r = 1.0D+00 do k = 1, nm r = r * ( c - a + k - 1.0D+00 ) * ( c - b + k - 1.0D+00 ) & / ( k * ( c + k - 1.0D+00 ) ) * x hf = hf + r end do hf = ( 1.0D+00 - x )**( c - a - b ) * hf return end if aa = a bb = b x1 = x if ( x < 0.0D+00 ) then x = x / ( x - 1.0D+00 ) if ( a < c .and. b < a .and. 0.0D+00 < b ) then a = bb b = aa end if b = c - b end if if ( 0.75D+00 <= x ) then gm = 0.0D+00 if ( abs ( c - a - b - aint ( c - a - b ) ) < 1.0D-15 ) then m = int ( c - a - b ) ga = gamma ( a ) gb = gamma ( b ) gc = gamma ( c ) gam = gamma ( a + m ) gbm = gamma ( b + m ) pa = r8_psi ( a ) pb = r8_psi ( b ) if ( m /= 0 ) then gm = 1.0D+00 end if do j = 1, abs ( m ) - 1 gm = gm * j end do rm = 1.0D+00 do j = 1, abs ( m ) rm = rm * j end do f0 = 1.0D+00 r0 = 1.0D+00 r1 = 1.0D+00 sp0 = 0.0D+00 sp = 0.0D+00 if ( 0 <= m ) then c0 = gm * gc / ( gam * gbm ) c1 = - gc * ( x - 1.0D+00 )**m / ( ga * gb * rm ) do k = 1, m - 1 r0 = r0 * ( a + k - 1.0D+00 ) * ( b + k - 1.0D+00 ) & / ( k * ( k - m ) ) * ( 1.0D+00 - x ) f0 = f0 + r0 end do do k = 1, m sp0 = sp0 + 1.0D+00 / ( a + k - 1.0D+00 ) & + 1.0D+00 / ( b + k - 1.0D+00 ) - 1.0D+00 / real ( k, kind = rk ) end do f1 = pa + pb + sp0 + 2.0D+00 * el + log ( 1.0D+00 - x ) hw = f1 do k = 1, 250 sp = sp + ( 1.0D+00 - a ) / ( k * ( a + k - 1.0D+00 ) ) & + ( 1.0D+00 - b ) / ( k * ( b + k - 1.0D+00 ) ) sm = 0.0D+00 do j = 1, m sm = sm + ( 1.0D+00 - a ) & / ( ( j + k ) * ( a + j + k - 1.0D+00 ) ) & + 1.0D+00 / ( b + j + k - 1.0D+00 ) end do rp = pa + pb + 2.0D+00 * el + sp + sm + log ( 1.0D+00 - x ) r1 = r1 * ( a + m + k - 1.0D+00 ) * ( b + m + k - 1.0D+00 ) & / ( k * ( m + k ) ) * ( 1.0D+00 - x ) f1 = f1 + r1 * rp if ( abs ( f1 - hw ) < abs ( f1 ) * eps ) then exit end if hw = f1 end do hf = f0 * c0 + f1 * c1 else if ( m < 0 ) then m = - m c0 = gm * gc / ( ga * gb * ( 1.0D+00 - x )**m ) c1 = - ( - 1 )**m * gc / ( gam * gbm * rm ) do k = 1, m - 1 r0 = r0 * ( a - m + k - 1.0D+00 ) * ( b - m + k - 1.0D+00 ) & / ( k * ( k - m ) ) * ( 1.0D+00 - x ) f0 = f0 + r0 end do do k = 1, m sp0 = sp0 + 1.0D+00 / real ( k, kind = rk ) end do f1 = pa + pb - sp0 + 2.0D+00 * el + log ( 1.0D+00 - x ) hw = f1 do k = 1, 250 sp = sp + ( 1.0D+00 - a ) & / ( k * ( a + k - 1.0D+00 ) ) & + ( 1.0D+00 - b ) / ( k * ( b + k - 1.0D+00 ) ) sm = 0.0D+00 do j = 1, m sm = sm + 1.0D+00 / real ( j + k, kind = rk ) end do rp = pa + pb + 2.0D+00 * el + sp - sm + log ( 1.0D+00 - x ) r1 = r1 * ( a + k - 1.0D+00 ) * ( b + k - 1.0D+00 ) & / ( k * ( m + k ) ) * ( 1.0D+00 - x ) f1 = f1 + r1 * rp if ( abs ( f1 - hw ) < abs ( f1 ) * eps ) then exit end if hw = f1 end do hf = f0 * c0 + f1 * c1 end if else ga = gamma ( a ) gb = gamma ( b ) gc = gamma ( c ) gca = gamma ( c - a ) gcb = gamma ( c - b ) gcab = gamma ( c - a - b ) gabc = gamma ( a + b - c ) c0 = gc * gcab / ( gca * gcb ) c1 = gc * gabc / ( ga * gb ) * ( 1.0D+00 - x )**( c - a - b ) hf = 0.0D+00 hw = hf r0 = c0 r1 = c1 do k = 1, 250 r0 = r0 * ( a + k - 1.0D+00 ) * ( b + k - 1.0D+00 ) & / ( k * ( a + b - c + k ) ) * ( 1.0D+00 - x ) r1 = r1 * ( c - a + k - 1.0D+00 ) * ( c - b + k - 1.0D+00 ) & / ( k * ( c - a - b + k ) ) * ( 1.0D+00 - x ) hf = hf + r0 + r1 if ( abs ( hf - hw ) < abs ( hf ) * eps ) then exit end if hw = hf end do hf = hf + c0 + c1 end if else a0 = 1.0D+00 if ( a < c .and. c < 2.0D+00 * a .and. b < c .and. c < 2.0D+00 * b ) then a0 = ( 1.0D+00 - x )**( c - a - b ) a = c - a b = c - b end if hf = 1.0D+00 hw = hf r = 1.0D+00 do k = 1, 250 r = r * ( a + k - 1.0D+00 ) * ( b + k - 1.0D+00 ) & / ( k * ( c + k - 1.0D+00 ) ) * x hf = hf + r if ( abs ( hf - hw ) <= abs ( hf ) * eps ) then exit end if hw = hf end do hf = a0 * hf end if if ( x1 < 0.0D+00 ) then x = x1 c0 = 1.0D+00 / ( 1.0D+00 - x )**aa hf = c0 * hf end if a = aa b = bb if ( 120 < k ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8_HYPER_2F1 - Warning!' write ( *, '(a)' ) ' A large number of iterations were needed.' write ( *, '(a)' ) ' The accuracy of the results should be checked.' end if return end function r8_psi ( xx ) !*****************************************************************************80 ! !! R8_PSI evaluates the function Psi(X). ! ! Discussion: ! ! This routine evaluates the logarithmic derivative of the ! Gamma function, ! ! PSI(X) = d/dX ( GAMMA(X) ) / GAMMA(X) ! = d/dX LN ( GAMMA(X) ) ! ! for real X, where either ! ! - XMAX1 < X < - XMIN, and X is not a negative integer, ! ! or ! ! XMIN < X. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! Original FORTRAN77 version by William Cody. ! FORTRAN90 version by John Burkardt. ! ! Reference: ! ! William Cody, Anthony Strecok, Henry Thacher, ! Chebyshev Approximations for the Psi Function, ! Mathematics of Computation, ! Volume 27, Number 121, January 1973, pages 123-127. ! ! Parameters: ! ! Input, real ( kind = rk ) XX, the argument of the function. ! ! Output, real ( kind = rk ) R8_PSI, the value of the function. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) real ( kind = rk ) aug real ( kind = rk ) den integer i integer n integer nq real ( kind = rk ), dimension ( 9 ) :: p1 = (/ & 4.5104681245762934160D-03, & 5.4932855833000385356D+00, & 3.7646693175929276856D+02, & 7.9525490849151998065D+03, & 7.1451595818951933210D+04, & 3.0655976301987365674D+05, & 6.3606997788964458797D+05, & 5.8041312783537569993D+05, & 1.6585695029761022321D+05 /) real ( kind = rk ), dimension ( 7 ) :: p2 = (/ & -2.7103228277757834192D+00, & -1.5166271776896121383D+01, & -1.9784554148719218667D+01, & -8.8100958828312219821D+00, & -1.4479614616899842986D+00, & -7.3689600332394549911D-02, & -6.5135387732718171306D-21 /) real ( kind = rk ), parameter :: piov4 = 0.78539816339744830962D+00 real ( kind = rk ), dimension ( 8 ) :: q1 = (/ & 9.6141654774222358525D+01, & 2.6287715790581193330D+03, & 2.9862497022250277920D+04, & 1.6206566091533671639D+05, & 4.3487880712768329037D+05, & 5.4256384537269993733D+05, & 2.4242185002017985252D+05, & 6.4155223783576225996D-08 /) real ( kind = rk ), dimension ( 6 ) :: q2 = (/ & 4.4992760373789365846D+01, & 2.0240955312679931159D+02, & 2.4736979003315290057D+02, & 1.0742543875702278326D+02, & 1.7463965060678569906D+01, & 8.8427520398873480342D-01 /) real ( kind = rk ) r8_psi real ( kind = rk ) sgn real ( kind = rk ) upper real ( kind = rk ) w real ( kind = rk ) x real ( kind = rk ), parameter :: x01 = 187.0D+00 real ( kind = rk ), parameter :: x01d = 128.0D+00 real ( kind = rk ), parameter :: x02 = 6.9464496836234126266D-04 real ( kind = rk ), parameter :: xinf = 1.70D+38 real ( kind = rk ), parameter :: xlarge = 2.04D+15 real ( kind = rk ), parameter :: xmax1 = 3.60D+16 real ( kind = rk ), parameter :: xmin1 = 5.89D-39 real ( kind = rk ), parameter :: xsmall = 2.05D-09 real ( kind = rk ) xx real ( kind = rk ) z x = xx w = abs ( x ) aug = 0.0D+00 ! ! Check for valid arguments, then branch to appropriate algorithm. ! if ( xmax1 <= - x .or. w < xmin1 ) then if ( 0.0D+00 < x ) then r8_psi = - xinf else r8_psi = xinf end if return end if if ( x < 0.5D+00 ) then ! ! X < 0.5, use reflection formula: psi(1-x) = psi(x) + pi * cot(pi*x) ! Use 1/X for PI*COTAN(PI*X) when XMIN1 < |X| <= XSMALL. ! if ( w <= xsmall ) then aug = - 1.0D+00 / x ! ! Argument reduction for cotangent. ! else if ( x < 0.0D+00 ) then sgn = piov4 else sgn = - piov4 end if w = w - real ( int ( w ), kind = rk ) nq = int ( w * 4.0D+00 ) w = 4.0D+00 * ( w - real ( nq, kind = rk ) * 0.25D+00 ) ! ! W is now related to the fractional part of 4.0 * X. ! Adjust argument to correspond to values in the first ! quadrant and determine the sign. ! n = nq / 2 if ( n + n /= nq ) then w = 1.0D+00 - w end if z = piov4 * w if ( mod ( n, 2 ) /= 0 ) then sgn = - sgn end if ! ! Determine the final value for -pi * cotan(pi*x). ! n = ( nq + 1 ) / 2 if ( mod ( n, 2 ) == 0 ) then ! ! Check for singularity. ! if ( z == 0.0D+00 ) then if ( 0.0D+00 < x ) then r8_psi = - xinf else r8_psi = xinf end if return end if aug = sgn * ( 4.0D+00 / tan ( z ) ) else aug = sgn * ( 4.0D+00 * tan ( z ) ) end if end if x = 1.0D+00 - x end if ! ! 0.5 <= X <= 3.0. ! if ( x <= 3.0D+00 ) then den = x upper = p1(1) * x do i = 1, 7 den = ( den + q1(i) ) * x upper = ( upper + p1(i+1) ) * x end do den = ( upper + p1(9) ) / ( den + q1(8) ) x = ( x - x01 / x01d ) - x02 r8_psi = den * x + aug return end if ! ! 3.0 < X. ! if ( x < xlarge ) then w = 1.0D+00 / ( x * x ) den = w upper = p2(1) * w do i = 1, 5 den = ( den + q2(i) ) * w upper = ( upper + p2(i+1) ) * w end do aug = ( upper + p2(7) ) / ( den + q2(6) ) - 0.5D+00 / x + aug end if r8_psi = aug + log ( x ) return end subroutine r8mat_data_read ( input_filename, m, n, table ) !*****************************************************************************80 ! !! R8MAT_DATA_READ reads data from an R8MAT file. ! ! Discussion: ! ! The file may contain more than N points, but this routine will ! return after reading N of them. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 18 October 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILENAME, the name of the input file. ! ! Input, integer M, the spatial dimension. ! ! Input, integer N, the number of points. ! ! Output, real ( kind = rk ) TABLE(M,N), the table data. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) integer m integer n integer ierror character ( len = * ) input_filename integer input_status integer input_unit integer j character ( len = 255 ) line real ( kind = rk ) table(m,n) real ( kind = rk ) x(m) ierror = 0 call get_unit ( input_unit ) open ( unit = input_unit, file = input_filename, status = 'old', & iostat = input_status ) if ( input_status /= 0 ) then ierror = 1 write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8MAT_DATA_READ - Fatal error!' write ( *, '(a,i8)' ) ' Could not open the input file "' // & trim ( input_filename ) // '" on unit ', input_unit stop end if j = 0 do while ( j < n ) read ( input_unit, '(a)', iostat = input_status ) line if ( input_status /= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8MAT_DATA_READ - Fatal error!' write ( *, '(a)' ) ' Error while reading lines of data.' write ( *, '(a,i8)' ) ' Number of values expected per line M = ', m write ( *, '(a,i8)' ) ' Number of data lines read, J = ', j write ( *, '(a,i8)' ) ' Number of data lines needed, N = ', n stop end if if ( line(1:1) == '#' .or. len_trim ( line ) == 0 ) then cycle end if call s_to_r8vec ( line, m, x, ierror ) if ( ierror /= 0 ) then cycle end if j = j + 1 table(1:m,j) = x(1:m) end do close ( unit = input_unit ) return end subroutine r8mat_header_read ( input_filename, m, n ) !*****************************************************************************80 ! !! R8MAT_HEADER_READ reads the header from an R8MAT file. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 07 September 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILENAME, the name of the input file. ! ! Output, integer M, spatial dimension. ! ! Output, integer N, the number of points. ! implicit none character ( len = * ) input_filename integer m integer n call file_column_count ( input_filename, m ) if ( m <= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8MAT_HEADER_READ - Fatal error!' write ( *, '(a)' ) ' There was some kind of I/O problem while trying' write ( *, '(a)' ) ' to count the number of data columns in' write ( *, '(a)' ) ' the file "' // trim ( input_filename ) // '".' stop end if call file_row_count ( input_filename, n ) if ( n <= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8MAT_HEADER_READ - Fatal error!' write ( *, '(a)' ) ' There was some kind of I/O problem while trying' write ( *, '(a)' ) ' to count the number of data rows in' write ( *, '(a)' ) ' the file "' // trim ( input_filename ) // '".' stop end if 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 MIT license. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) S, a string to be examined. ! ! Output, integer IVAL, the integer value read from the string. ! If the string is blank, then IVAL will be returned 0. ! ! Output, integer IERROR, an error flag. ! 0, no error. ! 1, an error occurred. ! ! Output, integer LENGTH, the number of characters of S ! used to make IVAL. ! implicit none character c integer i integer ierror integer isgn integer istate integer ival integer 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 MIT license. ! ! Modified: ! ! 11 February 2008 ! ! 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 = rk ) DVAL, the value read from the string. ! ! Output, integer 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 LENGTH, the number of characters read ! to form the number, including any terminating ! characters such as a trailing comma or blanks. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) logical ch_eqi character c real ( kind = rk ) dval integer ierror integer ihave integer isgn integer iterm integer jbot integer jsgn integer jtop integer length integer nchar integer ndig real ( kind = rk ) rbot real ( kind = rk ) rexp real ( kind = rk ) 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 = rk ) else if ( ihave == 5 ) then rtop = 10.0D+00 * rtop + real ( ndig, kind = rk ) 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 = rk ) & / real ( jbot, kind = rk ) ) end if end if dval = real ( isgn, kind = rk ) * rexp * rtop / rbot return end subroutine s_to_r8vec ( s, n, rvec, ierror ) !*****************************************************************************80 ! !! S_TO_R8VEC reads an R8VEC from a string. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) S, the string to be read. ! ! Input, integer N, the number of values expected. ! ! Output, real ( kind = rk ) RVEC(N), the values read from the string. ! ! Output, integer IERROR, error flag. ! 0, no errors occurred. ! -K, could not read data for entries -K through N. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) integer n integer i integer ierror integer ilo integer lchar real ( kind = rk ) rvec(n) character ( len = * ) s i = 0 ierror = 0 ilo = 1 do while ( i < n ) i = i + 1 call s_to_r8 ( s(ilo:), rvec(i), ierror, lchar ) if ( ierror /= 0 ) then ierror = -i exit end if ilo = ilo + lchar end do return end subroutine s_word_count ( s, nword ) !*****************************************************************************80 ! !! S_WORD_COUNT counts the number of "words" in a string. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 11 February 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) S, the string to be examined. ! ! Output, integer NWORD, the number of "words" in the string. ! Words are presumed to be separated by one or more blanks. ! implicit none logical blank integer i integer lens integer nword character ( len = * ) s nword = 0 lens = len ( s ) if ( lens <= 0 ) then return end if blank = .true. do i = 1, lens if ( s(i:i) == ' ' ) then blank = .true. else if ( blank ) then nword = nword + 1 blank = .false. end if end do 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: ! ! 18 May 2013 ! ! 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,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