subroutine burgers_viscous_time_exact1 ( nu, vxn, vx, vtn, vt, vu ) !*****************************************************************************80 ! !! BURGERS_VISCOUS_TIME_EXACT1 evaluates solution #1 to the Burgers equation. ! ! Discussion: ! ! The form of the Burgers equation considered here is ! ! du du d^2 u ! -- + u * -- = nu * ----- ! dt dx dx^2 ! ! for -1.0 < x < +1.0, and 0 < t. ! ! Initial conditions are u(x,0) = - sin(pi*x). Boundary conditions ! are u(-1,t) = u(+1,t) = 0. The viscosity parameter nu is taken ! to be 0.01 / pi, although this is not essential. ! ! The authors note an integral representation for the solution u(x,t), ! and present a better version of the formula that is amenable to ! approximation using Hermite quadrature. ! ! This program library does little more than evaluate the exact solution ! at a user-specified set of points, using the quadrature rule. ! Internally, the order of this quadrature rule is set to 8, but the ! user can easily modify this value if greater accuracy is desired. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 25 September 2015 ! ! Author: ! ! John Burkardt. ! ! Reference: ! ! Claude Basdevant, Michel Deville, Pierre Haldenwang, J Lacroix, ! J Ouazzani, Roger Peyret, Paolo Orlandi, Anthony Patera, ! Spectral and finite difference solutions of the Burgers equation, ! Computers and Fluids, ! Volume 14, Number 1, 1986, pages 23-41. ! ! Parameters: ! ! Input, real ( kind = 8 ) NU, the viscosity. ! ! Input, integer ( kind = 4 ) VXN, the number of spatial grid points. ! ! Input, real ( kind = 8 ) VX(VXN), the spatial grid points. ! ! Input, integer ( kind = 4 ) VTN, the number of time grid points. ! ! Input, real ( kind = 8 ) VT(VTN), the time grid points. ! ! Output, real ( kind = 8 ) VU(VXN,VTN), the solution of the Burgers ! equation at each space and time grid point. ! implicit none integer ( kind = 4 ), parameter :: qn = 8 integer ( kind = 4 ) vtn integer ( kind = 4 ) vxn real ( kind = 8 ) bot real ( kind = 8 ) c real ( kind = 8 ) nu integer ( kind = 4 ) qi real ( kind = 8 ) qw(qn) real ( kind = 8 ) qx(qn) real ( kind = 8 ), parameter :: r8_pi = 3.141592653589793D+00 real ( kind = 8 ) vt(vtn) integer ( kind = 4 ) vti real ( kind = 8 ) vx(vxn) integer ( kind = 4 ) vxi real ( kind = 8 ) vu(vxn,vtn) real ( kind = 8 ) top ! ! Compute the rule. ! call hermite_ek_compute ( qn, qx, qw ) ! ! Evaluate U(X,T) for later times. ! do vti = 1, vtn if ( vt(vti) == 0.0D+00 ) then do vxi = 1, vxn vu(vxi,vti) = - sin ( r8_pi * vx(vxi) ) end do else do vxi = 1, vxn top = 0.0D+00 bot = 0.0D+00 do qi = 1, qn c = 2.0D+00 * sqrt ( nu * vt(vti) ) top = top - qw(qi) * c * sin ( r8_pi * ( vx(vxi) - c * qx(qi) ) ) & * exp ( - cos ( r8_pi * ( vx(vxi) - c * qx(qi) ) ) & / ( 2.0D+00 * r8_pi * nu ) ) bot = bot + qw(qi) * c & * exp ( - cos ( r8_pi * ( vx(vxi) - c * qx(qi) ) ) & / ( 2.0D+00 * r8_pi * nu ) ) vu(vxi,vti) = top / bot end do end do end if end do return end subroutine burgers_viscous_time_exact2 ( nu, xn, x, tn, t, u ) !*****************************************************************************80 ! !! BURGERS_VISCOUS_TIME_EXACT2 evaluates solution #2 to the Burgers equation. ! ! Discussion: ! ! The form of the Burgers equation considered here is ! ! du du d^2 u ! -- + u * -- = nu * ----- ! dt dx dx^2 ! ! for 0.0 < x < 2 Pi and 0 < t. ! ! The initial condition is ! ! u(x,0) = 4 - 2 * nu * dphi(x,0)/dx / phi(x,0) ! ! where ! ! phi(x,t) = exp ( - ( x-4*t ) / ( 4*nu*(t+1) ) ) ! + exp ( - ( x-4*t-2*pi ) / ( 4*nu*(t+1) ) ) ! ! The boundary conditions are periodic: ! ! u(0,t) = u(2 Pi,t) ! ! The viscosity parameter nu may be taken to be 0.01, but other values ! may be chosen. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 25 September 2015 ! ! Author: ! ! John Burkardt. ! ! Parameters: ! ! Input, real ( kind = 8 ) NU, the viscosity. ! ! Input, integer ( kind = 4 ) XN, the number of spatial grid points. ! ! Input, real ( kind = 8 ) X(XN), the spatial grid points. ! ! Input, integer ( kind = 4 ) TN, the number of time grid points. ! ! Input, real ( kind = 8 ) T(TN), the time grid points. ! ! Output, real ( kind = 8 ) U(XN,TN), the solution of the Burgers ! equation at each space and time grid point. ! implicit none integer ( kind = 4 ) tn integer ( kind = 4 ) xn real ( kind = 8 ) a real ( kind = 8 ) b real ( kind = 8 ) c real ( kind = 8 ) dphi integer ( kind = 4 ) i integer ( kind = 4 ) j real ( kind = 8 ) nu real ( kind = 8 ) phi real ( kind = 8 ), parameter :: r8_pi = 3.141592653589793D+00 real ( kind = 8 ) t(tn) real ( kind = 8 ) u(xn,tn) real ( kind = 8 ) x(xn) do j = 1, tn do i = 1, xn a = ( x(i) - 4.0D+00 * t(j) ) b = ( x(i) - 4.0D+00 * t(j) - 2.0D+00 * r8_pi ) c = 4.0D+00 * nu * ( t(j) + 1.0D+00 ) phi = exp ( - a ** 2 / c ) + exp ( - b ** 2 / c ) dphi = - 2.0D+00 * a * exp ( - a ** 2 / c ) / c & - 2.0D+00 * b * exp ( - b ** 2 / c ) / c u(i,j) = 4.0D+00 - 2.0D+00 * nu * dphi / phi end do end do return end subroutine get_unit ( iunit ) !*****************************************************************************80 ! !! GET_UNIT returns a free FORTRAN unit number. ! ! Discussion: ! ! A "free" FORTRAN unit number is a value 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 a value 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: ! ! 26 October 2008 ! ! 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 hermite_ek_compute ( n, x, w ) !*****************************************************************************80 ! !! HERMITE_EK_COMPUTE computes a Gauss-Hermite quadrature rule. ! ! Discussion: ! ! The code uses an algorithm by Elhay and Kautsky. ! ! The abscissas are the zeros of the N-th order Hermite polynomial. ! ! The integral: ! ! integral ( -oo < x < +oo ) exp ( - x * x ) * f(x) dx ! ! The quadrature rule: ! ! sum ( 1 <= i <= n ) w(i) * f ( x(i) ) ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 15 April 2011 ! ! Author: ! ! Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky. ! FORTRAN90 version by John Burkardt. ! ! Reference: ! ! Sylvan Elhay, Jaroslav Kautsky, ! Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of ! Interpolatory Quadrature, ! ACM Transactions on Mathematical Software, ! Volume 13, Number 4, December 1987, pages 399-415. ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, the number of abscissas. ! ! Output, real ( kind = 8 ) X(N), the abscissas. ! ! Output, real ( kind = 8 ) W(N), the weights. ! implicit none integer ( kind = 4 ) n real ( kind = 8 ) bj(n) integer ( kind = 4 ) i real ( kind = 8 ) r8_gamma real ( kind = 8 ) w(n) real ( kind = 8 ) x(n) real ( kind = 8 ) zemu ! ! Define the zero-th moment. ! zemu = r8_gamma ( 1.0D+00 / 2.0D+00 ) ! ! Define the Jacobi matrix. ! do i = 1, n bj(i) = real ( i, kind = 8 ) / 2.0D+00 end do bj(1:n) = sqrt ( bj(1:n) ) x(1:n) = 0.0D+00 w(1) = sqrt ( zemu ) w(2:n) = 0.0D+00 ! ! Diagonalize the Jacobi matrix. ! call imtqlx ( n, x, bj, w ) w(1:n) = w(1:n)**2 return end subroutine imtqlx ( n, d, e, z ) !*****************************************************************************80 ! !! IMTQLX diagonalizes a symmetric tridiagonal matrix. ! ! Discussion: ! ! This routine is a slightly modified version of the EISPACK routine to ! perform the implicit QL algorithm on a symmetric tridiagonal matrix. ! ! The authors thank the authors of EISPACK for permission to use this ! routine. ! ! It has been modified to produce the product Q' * Z, where Z is an input ! vector and Q is the orthogonal matrix diagonalizing the input matrix. ! The changes consist (essentially) of applying the orthogonal ! transformations directly to Z as they are generated. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 27 December 2009 ! ! Author: ! ! Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky. ! FORTRAN90 version by John Burkardt. ! ! Reference: ! ! Sylvan Elhay, Jaroslav Kautsky, ! Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of ! Interpolatory Quadrature, ! ACM Transactions on Mathematical Software, ! Volume 13, Number 4, December 1987, pages 399-415. ! ! Roger Martin, James Wilkinson, ! The Implicit QL Algorithm, ! Numerische Mathematik, ! Volume 12, Number 5, December 1968, pages 377-383. ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, the order of the matrix. ! ! Input/output, real ( kind = 8 ) D(N), the diagonal entries of the matrix. ! On output, the information in D has been overwritten. ! ! Input/output, real ( kind = 8 ) E(N), the subdiagonal entries of the ! matrix, in entries E(1) through E(N-1). On output, the information in ! E has been overwritten. ! ! Input/output, real ( kind = 8 ) Z(N). On input, a vector. On output, ! the value of Q' * Z, where Q is the matrix that diagonalizes the ! input symmetric tridiagonal matrix. ! implicit none integer ( kind = 4 ) n real ( kind = 8 ) b real ( kind = 8 ) c real ( kind = 8 ) d(n) real ( kind = 8 ) e(n) real ( kind = 8 ) f real ( kind = 8 ) g integer ( kind = 4 ) i integer ( kind = 4 ) ii integer ( kind = 4 ), parameter :: itn = 30 integer ( kind = 4 ) j integer ( kind = 4 ) k integer ( kind = 4 ) l integer ( kind = 4 ) m integer ( kind = 4 ) mml real ( kind = 8 ) p real ( kind = 8 ) prec real ( kind = 8 ) r real ( kind = 8 ) s real ( kind = 8 ) z(n) prec = epsilon ( prec ) if ( n == 1 ) then return end if e(n) = 0.0D+00 do l = 1, n j = 0 do do m = l, n if ( m == n ) then exit end if if ( abs ( e(m) ) <= prec * ( abs ( d(m) ) + abs ( d(m+1) ) ) ) then exit end if end do p = d(l) if ( m == l ) then exit end if if ( itn <= j ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'IMTQLX - Fatal error!' write ( *, '(a)' ) ' Iteration limit exceeded.' write ( *, '(a,i8)' ) ' J = ', j write ( *, '(a,i8)' ) ' L = ', l write ( *, '(a,i8)' ) ' M = ', m write ( *, '(a,i8)' ) ' N = ', n stop end if j = j + 1 g = ( d(l+1) - p ) / ( 2.0D+00 * e(l) ) r = sqrt ( g * g + 1.0D+00 ) g = d(m) - p + e(l) / ( g + sign ( r, g ) ) s = 1.0D+00 c = 1.0D+00 p = 0.0D+00 mml = m - l do ii = 1, mml i = m - ii f = s * e(i) b = c * e(i) if ( abs ( g ) <= abs ( f ) ) then c = g / f r = sqrt ( c * c + 1.0D+00 ) e(i+1) = f * r s = 1.0D+00 / r c = c * s else s = f / g r = sqrt ( s * s + 1.0D+00 ) e(i+1) = g * r c = 1.0D+00 / r s = s * c end if g = d(i+1) - p r = ( d(i) - g ) * s + 2.0D+00 * c * b p = s * r d(i+1) = g + p g = c * r - b f = z(i+1) z(i+1) = s * z(i) + c * f z(i) = c * z(i) - s * f end do d(l) = d(l) - p e(l) = g e(m) = 0.0D+00 end do end do ! ! Sorting. ! do ii = 2, n i = ii - 1 k = i p = d(i) do j = ii, n if ( d(j) < p ) then k = j p = d(j) end if end do if ( k /= i ) then d(k) = d(i) d(i) = p p = z(i) z(i) = z(k) z(k) = p end if end do return end function r8_gamma ( x ) !*****************************************************************************80 ! !! R8_GAMMA evaluates Gamma(X) for a real argument. ! ! Discussion: ! ! 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: ! ! Original FORTRAN77 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 integer ( kind = 4 ) i integer ( kind = 4 ) n 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 ) 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 parity = .false. fact = 1.0D+00 n = 0 y = x ! ! Argument is negative. ! if ( y <= 0.0D+00 ) then y = - x y1 = aint ( y ) res = y - y1 if ( res /= 0.0D+00 ) then if ( y1 /= aint ( y1 * 0.5D+00 ) * 2.0D+00 ) then parity = .true. end if fact = - pi / sin ( pi * res ) y = y + 1.0D+00 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 = 1.0D+00 / y else res = xinf r8_gamma = res return end if else if ( y < 12.0D+00 ) then y1 = y ! ! 0.0 < argument < 1.0. ! if ( y < 1.0D+00 ) then z = y y = y + 1.0D+00 ! ! 1.0 < argument < 12.0. ! Reduce argument if necessary. ! else n = int ( y ) - 1 y = y - real ( n, kind = 8 ) z = y - 1.0D+00 end if ! ! Evaluate approximation for 1.0 < argument < 2.0. ! xnum = 0.0D+00 xden = 1.0D+00 do i = 1, 8 xnum = ( xnum + p(i) ) * z xden = xden * z + q(i) end do res = xnum / xden + 1.0D+00 ! ! 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 + 1.0D+00 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 - 0.5D+00 ) * 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 /= 1.0D+00 ) then res = fact / res end if r8_gamma = res return end subroutine r8mat_print ( m, n, a, title ) !*****************************************************************************80 ! !! R8MAT_PRINT prints an R8MAT. ! ! Discussion: ! ! An R8MAT is an array of R8 values. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 12 September 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) M, the number of rows in A. ! ! Input, integer ( kind = 4 ) N, the number of columns in A. ! ! Input, real ( kind = 8 ) A(M,N), the matrix. ! ! Input, character ( len = * ) TITLE, a title. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n real ( kind = 8 ) a(m,n) character ( len = * ) title call r8mat_print_some ( m, n, a, 1, 1, m, n, title ) return end subroutine r8mat_print_some ( m, n, a, ilo, jlo, ihi, jhi, title ) !*****************************************************************************80 ! !! R8MAT_PRINT_SOME prints some of an R8MAT. ! ! Discussion: ! ! An R8MAT is an array of R8 values. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 10 September 2009 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) M, N, the number of rows and columns. ! ! Input, real ( kind = 8 ) A(M,N), an M by N matrix to be printed. ! ! Input, integer ( kind = 4 ) ILO, JLO, the first row and column to print. ! ! Input, integer ( kind = 4 ) IHI, JHI, the last row and column to print. ! ! Input, character ( len = * ) TITLE, a title. ! implicit none integer ( kind = 4 ), parameter :: incx = 5 integer ( kind = 4 ) m integer ( kind = 4 ) n real ( kind = 8 ) a(m,n) character ( len = 14 ) ctemp(incx) integer ( kind = 4 ) i integer ( kind = 4 ) i2hi integer ( kind = 4 ) i2lo integer ( kind = 4 ) ihi integer ( kind = 4 ) ilo integer ( kind = 4 ) inc integer ( kind = 4 ) j integer ( kind = 4 ) j2 integer ( kind = 4 ) j2hi integer ( kind = 4 ) j2lo integer ( kind = 4 ) jhi integer ( kind = 4 ) jlo character ( len = * ) title write ( *, '(a)' ) ' ' write ( *, '(a)' ) trim ( title ) if ( m <= 0 .or. n <= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' (None)' return end if do j2lo = max ( jlo, 1 ), min ( jhi, n ), incx j2hi = j2lo + incx - 1 j2hi = min ( j2hi, n ) j2hi = min ( j2hi, jhi ) inc = j2hi + 1 - j2lo write ( *, '(a)' ) ' ' do j = j2lo, j2hi j2 = j + 1 - j2lo write ( ctemp(j2), '(i8,6x)' ) j end do write ( *, '('' Col '',5a14)' ) ctemp(1:inc) write ( *, '(a)' ) ' Row' write ( *, '(a)' ) ' ' i2lo = max ( ilo, 1 ) i2hi = min ( ihi, m ) do i = i2lo, i2hi do j2 = 1, inc j = j2lo - 1 + j2 if ( a(i,j) == real ( int ( a(i,j) ), kind = 8 ) ) then write ( ctemp(j2), '(f8.0,6x)' ) a(i,j) else write ( ctemp(j2), '(g14.6)' ) a(i,j) end if end do write ( *, '(i5,a,5a14)' ) i, ':', ( ctemp(j), j = 1, inc ) end do end do return end subroutine r8mat_write ( output_filename, m, n, table ) !*****************************************************************************80 ! !! R8MAT_WRITE writes an R8MAT file. ! ! Discussion: ! ! An R8MAT is an array of R8 values. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 31 May 2009 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) OUTPUT_FILENAME, 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 data. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n integer ( kind = 4 ) j character ( len = * ) output_filename integer ( kind = 4 ) output_status integer ( kind = 4 ) output_unit character ( len = 30 ) string real ( kind = 8 ) table(m,n) ! ! Open the file. ! call get_unit ( output_unit ) open ( unit = output_unit, file = output_filename, & status = 'replace', iostat = output_status ) if ( output_status /= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8MAT_WRITE - Fatal error!' write ( *, '(a,i8)' ) ' Could not open the output file "' // & trim ( output_filename ) // '" on unit ', output_unit output_unit = -1 stop end if ! ! Create a format string. ! ! For less precision in the output file, try: ! ! '(', m, 'g', 14, '.', 6, ')' ! if ( 0 < m .and. 0 < n ) then write ( string, '(a1,i8,a1,i8,a1,i8,a1)' ) '(', m, 'g', 24, '.', 16, ')' ! ! Write the data. ! do j = 1, n write ( output_unit, string ) table(1:m,j) end do end if ! ! Close the file. ! close ( unit = output_unit ) return end subroutine r8vec_even ( n, alo, ahi, a ) !*****************************************************************************80 ! !! R8VEC_EVEN returns an R8VEC of evenly spaced values. ! ! Discussion: ! ! An R8VEC is a vector of R8's. ! ! If N is 1, then the midpoint is returned. ! ! Otherwise, the two endpoints are returned, and N-2 evenly ! spaced points between them. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 09 December 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, the number of values. ! ! Input, real ( kind = 8 ) ALO, AHI, the low and high values. ! ! Output, real ( kind = 8 ) A(N), N evenly spaced values. ! Normally, A(1) = ALO and A(N) = AHI. ! However, if N = 1, then A(1) = 0.5*(ALO+AHI). ! implicit none integer ( kind = 4 ) n real ( kind = 8 ) a(n) real ( kind = 8 ) ahi real ( kind = 8 ) alo integer ( kind = 4 ) i if ( n == 1 ) then a(1) = 0.5D+00 * ( alo + ahi ) else do i = 1, n a(i) = ( real ( n - i, kind = 8 ) * alo & + real ( i - 1, kind = 8 ) * ahi ) & / real ( n - 1, kind = 8 ) end do end if return end subroutine r8vec_print ( n, a, title ) !*****************************************************************************80 ! !! R8VEC_PRINT prints an R8VEC. ! ! Discussion: ! ! An R8VEC is a vector of R8's. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 22 August 2000 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, the number of components of the vector. ! ! Input, real ( kind = 8 ) A(N), the vector to be printed. ! ! Input, character ( len = * ) TITLE, a title. ! implicit none integer ( kind = 4 ) n real ( kind = 8 ) a(n) integer ( kind = 4 ) i character ( len = * ) title write ( *, '(a)' ) ' ' write ( *, '(a)' ) trim ( title ) write ( *, '(a)' ) ' ' do i = 1, n write ( *, '(2x,i8,a,1x,g16.8)' ) i, ':', a(i) 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 GNU LGPL license. ! ! Modified: ! ! 06 August 2005 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! None ! implicit none character ( len = 8 ) ampm integer ( kind = 4 ) d 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 integer ( kind = 4 ) values(8) integer ( kind = 4 ) 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