subroutine monomial_value ( m, n, e, x, v ) !*****************************************************************************80 ! !! MONOMIAL_VALUE evaluates a monomial. ! ! Discussion: ! ! F(X) = product ( 1 <= I <= M ) X(I)^EXPON(I) ! ! with the convention that 0^0 = 1. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 13 April 2009 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) M, the spatial dimension. ! ! Input, integer ( kind = 4 ) N, the number of points. ! ! Input, integer ( kind = 4 ) E(M), the exponents. ! ! Input, real ( kind = 8 ) X(M,N), the evaluation points. ! ! Output, real ( kind = 8 ) V(N), the monomial values. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n integer ( kind = 4 ) i integer ( kind = 4 ) e(m) real ( kind = 8 ) v(n) real ( kind = 8 ) x(m,n) v(1:n) = 1.0D+00 do i = 1, m if ( e(i) /= 0.0D+00 ) then v(1:n) = v(1:n) * x(i,1:n) ** e(i) end if end do return end subroutine pyramid01_monomial_integral ( expon, value ) !*****************************************************************************80 ! !! PYRAMID01_MONOMIAL_INTEGRAL: monomial integral in a unit pyramid. ! ! Discussion: ! ! This routine returns the integral of ! ! product ( 1 <= I <= 3 ) X(I)^EXPON(I) ! ! over the unit pyramid. ! ! The integration region is: ! ! - ( 1 - Z ) <= X <= 1 - Z ! - ( 1 - Z ) <= Y <= 1 - Z ! 0 <= Z <= 1. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 24 March 2008 ! ! Author: ! ! John Burkardt ! ! Reference: ! ! Arthur Stroud, ! Approximate Calculation of Multiple Integrals, ! Prentice Hall, 1971, ! ISBN: 0130438936, ! LC: QA311.S85. ! ! Parameters: ! ! Input, integer ( kind = 4 ) EXPON(3), the exponents. ! ! Output, real ( kind = 8 ) VALUE, the integral of the monomial. ! implicit none integer ( kind = 4 ) expon(3) integer ( kind = 4 ) i integer ( kind = 4 ) i_hi real ( kind = 8 ) r8_choose real ( kind = 8 ) r8_mop real ( kind = 8 ) value value = 0.0D+00 if ( mod ( expon(1), 2 ) == 0 .and. mod ( expon(2), 2 ) == 0 ) then i_hi = 2 + expon(1) + expon(2) do i = 0, i_hi value = value + r8_mop ( i ) * r8_choose ( i_hi, i ) & / real ( i + expon(3) + 1, kind = 8 ) end do value = value & * 2.0D+00 / real ( expon(1) + 1, kind = 8 ) & * 2.0D+00 / real ( expon(2) + 1, kind = 8 ) end if return end subroutine pyramid01_sample ( n, seed, x ) !*****************************************************************************80 ! !! PYRAMID01_SAMPLE: sample the unit pyramid. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 13 April 2014 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, the number of samples desired. ! ! Input/output, integer ( kind = 4 ) SEED, a seed for the random ! number generator. ! ! Output, real ( kind = 8 ) X(3,N), the sample values. ! implicit none integer ( kind = 4 ) n real ( kind = 8 ), parameter :: one_third = 1.0D+00 / 3.0D+00 integer ( kind = 4 ) seed real ( kind = 8 ) x(3,n) call r8mat_uniform_01 ( 3, n, seed, x ) x(3,1:n) = 1.0D+00 - x(3,1:n) ** one_third x(2,1:n) = ( 1.0D+00 - x(3,1:n) ) * ( 2.0D+00 * x(2,1:n) - 1.0D+00 ) x(1,1:n) = ( 1.0D+00 - x(3,1:n) ) * ( 2.0D+00 * x(1,1:n) - 1.0D+00 ) return end function pyramid01_volume ( ) !*****************************************************************************80 ! !! PYRAMID01_VOLUME: volume of a unit pyramid with square base. ! ! Discussion: ! ! The volume of this unit pyramid is 4/3. ! ! The integration region is: ! ! - ( 1 - Z ) <= X <= 1 - Z ! - ( 1 - Z ) <= Y <= 1 - Z ! 0 <= Z <= 1. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 22 March 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Output, real ( kind = 8 ) PYRAMID01_VOLUME, the volume. ! implicit none real ( kind = 8 ) pyramid01_volume pyramid01_volume = 4.0D+00 / 3.0D+00 return end function r8_choose ( n, k ) !*****************************************************************************80 ! !! R8_CHOOSE computes the binomial coefficient C(N,K) as an R8. ! ! Discussion: ! ! The value is calculated in such a way as to avoid overflow and ! roundoff. The calculation is done in R8 arithmetic. ! ! The formula used is: ! ! C(N,K) = N! / ( K! * (N-K)! ) ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 24 March 2008 ! ! Author: ! ! John Burkardt ! ! Reference: ! ! ML Wolfson, HV Wright, ! Algorithm 160: ! Combinatorial of M Things Taken N at a Time, ! Communications of the ACM, ! Volume 6, Number 4, April 1963, page 161. ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, K, are the values of N and K. ! ! Output, real ( kind = 8 ) R8_CHOOSE, the number of combinations of N ! things taken K at a time. ! implicit none integer ( kind = 4 ) i integer ( kind = 4 ) k integer ( kind = 4 ) mn integer ( kind = 4 ) mx integer ( kind = 4 ) n real ( kind = 8 ) r8_choose real ( kind = 8 ) value mn = min ( k, n - k ) if ( mn < 0 ) then value = 0.0D+00 else if ( mn == 0 ) then value = 1.0D+00 else mx = max ( k, n - k ) value = real ( mx + 1, kind = 8 ) do i = 2, mn value = ( value * real ( mx + i, kind = 8 ) ) / real ( i, kind = 8 ) end do end if r8_choose = value return end function r8_mop ( i ) !*****************************************************************************80 ! !! R8_MOP returns the I-th power of -1 as an R8. ! ! Discussion: ! ! An R8 is a real ( kind = 8 ) value. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 07 November 2007 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) I, the power of -1. ! ! Output, real ( kind = 8 ) R8_MOP, the I-th power of -1. ! implicit none integer ( kind = 4 ) i real ( kind = 8 ) r8_mop if ( mod ( i, 2 ) == 0 ) then r8_mop = + 1.0D+00 else r8_mop = - 1.0D+00 end if return end subroutine r8mat_transpose_print ( m, n, a, title ) !*****************************************************************************80 ! !! R8MAT_TRANSPOSE_PRINT prints an R8MAT, transposed. ! ! Discussion: ! ! An R8MAT is an MxN array of R8's, stored by (I,J) -> [I+J*M]. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 14 June 2004 ! ! 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, 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_transpose_print_some ( m, n, a, 1, 1, m, n, title ) return end subroutine r8mat_transpose_print_some ( m, n, a, ilo, jlo, ihi, jhi, title ) !*****************************************************************************80 ! !! R8MAT_TRANSPOSE_PRINT_SOME prints some of an R8MAT, transposed. ! ! Discussion: ! ! An R8MAT is an MxN array of R8's, stored by (I,J) -> [I+J*M]. ! ! 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 ) i2 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 ) 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 i2lo = max ( ilo, 1 ), min ( ihi, m ), incx i2hi = i2lo + incx - 1 i2hi = min ( i2hi, m ) i2hi = min ( i2hi, ihi ) inc = i2hi + 1 - i2lo write ( *, '(a)' ) ' ' do i = i2lo, i2hi i2 = i + 1 - i2lo write ( ctemp(i2), '(i8,6x)' ) i end do write ( *, '('' Row '',5a14)' ) ctemp(1:inc) write ( *, '(a)' ) ' Col' write ( *, '(a)' ) ' ' j2lo = max ( jlo, 1 ) j2hi = min ( jhi, n ) do j = j2lo, j2hi do i2 = 1, inc i = i2lo - 1 + i2 write ( ctemp(i2), '(g14.6)' ) a(i,j) end do write ( *, '(i5,a,5a14)' ) j, ':', ( ctemp(i), i = 1, inc ) end do end do return end subroutine r8mat_uniform_01 ( m, n, seed, r ) !*****************************************************************************80 ! !! R8MAT_UNIFORM_01 fills an R8MAT with unit pseudorandom numbers. ! ! Discussion: ! ! An R8MAT is an MxN array of R8's, stored by (I,J) -> [I+J*M]. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 11 August 2004 ! ! Author: ! ! John Burkardt ! ! Reference: ! ! Paul Bratley, Bennett Fox, Linus Schrage, ! A Guide to Simulation, ! Springer Verlag, pages 201-202, 1983. ! ! Bennett Fox, ! Algorithm 647: ! Implementation and Relative Efficiency of Quasirandom ! Sequence Generators, ! ACM Transactions on Mathematical Software, ! Volume 12, Number 4, pages 362-376, 1986. ! ! Peter Lewis, Allen Goodman, James Miller, ! A Pseudo-Random Number Generator for the System/360, ! IBM Systems Journal, ! Volume 8, pages 136-143, 1969. ! ! Parameters: ! ! Input, integer ( kind = 4 ) M, N, the number of rows and columns in ! the array. ! ! Input/output, integer ( kind = 4 ) SEED, the "seed" value, which ! should NOT be 0. On output, SEED has been updated. ! ! Output, real ( kind = 8 ) R(M,N), the array of pseudorandom values. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n integer ( kind = 4 ) i integer ( kind = 4 ), parameter :: i4_huge = 2147483647 integer ( kind = 4 ) j integer ( kind = 4 ) k integer ( kind = 4 ) seed real ( kind = 8 ) r(m,n) do j = 1, n do i = 1, m k = seed / 127773 seed = 16807 * ( seed - k * 127773 ) - k * 2836 if ( seed < 0 ) then seed = seed + i4_huge end if r(i,j) = real ( seed, kind = 8 ) * 4.656612875D-10 end do 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: ! ! 18 May 2013 ! ! 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