subroutine brownian_displacement_display ( k, n, m, d, t, dsq, header ) !*****************************************************************************80 ! !! BROWNIAN_DISPLACEMENT_DISPLAY displays average Brownian motion displacement. ! ! Discussion: ! ! Thanks to Feifei Xu for pointing out a missing factor of 2 in the ! displacement calculation. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 08 September 2018 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) K, the number of repetitions. ! ! Input, integer ( kind = 4 ) N, the number of time steps. ! ! Input, integer ( kind = 4 ) M, the spatial dimension. ! ! Input, real ( kind = 8 ) D, the diffusion coefficient. ! ! Input, real ( kind = 8 ) T, the total time. ! ! Input, real ( kind = 8 ) DSQ(K,N), the displacements over time for ! each repetition. ! ! Input, character ( len = * ) HEADER, an identifier for the output files. ! implicit none integer ( kind = 4 ) k integer ( kind = 4 ) n character ( len = 255 ) command_filename integer ( kind = 4 ) command_unit real ( kind = 8 ) d character ( len = 255 ) data_filename integer ( kind = 4 ) data_unit real ( kind = 8 ) dsq(k,n) real ( kind = 8 ) dsq_ave real ( kind = 8 ) dsq_ideal character ( len = * ) header integer ( kind = 4 ) i4_uniform integer ( kind = 4 ) ii(5) integer ( kind = 4 ) j integer ( kind = 4 ) m integer ( kind = 4 ) seed real ( kind = 8 ) t real ( kind = 8 ) ti seed = 123456789 ! ! Choose 5 paths at random. ! do j = 1, 5 ii(j) = i4_uniform ( 1, k, seed ) end do ! ! Create the data file. ! call get_unit ( data_unit ) data_filename = trim ( header ) // '_data.txt' open ( unit = data_unit, file = data_filename, status = 'replace' ) do j = 1, n ti = real ( j - 1, kind = 8 ) * t / real ( n - 1, kind = 8 ) dsq_ave = sum ( dsq(1:k,j) ) / real ( k, kind = 8 ) dsq_ideal = 2.0D+00 * real ( m, kind = 8 ) * d * ti write ( data_unit, '(8(2x,g14.6))' ) & ti, dsq(ii(1:5),j), dsq_ave, dsq_ideal end do close ( unit = data_unit ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' BROWNIAN_DISPLACEMENT data stored in "' & // trim ( data_filename ) // '".' ! ! Create the command file. ! call get_unit ( command_unit ) command_filename = trim ( header ) // '_commands.txt' open ( unit = command_unit, file = command_filename, status = 'replace' ) write ( command_unit, '(a)' ) '# ' // trim ( command_filename ) write ( command_unit, '(a)' ) '#' write ( command_unit, '(a)' ) '# Usage:' write ( command_unit, '(a)' ) '# gnuplot < ' // trim ( command_filename ) write ( command_unit, '(a)' ) '#' write ( command_unit, '(a)' ) 'set term png' write ( command_unit, '(a)' ) 'set output "' // trim ( header ) // '.png"' write ( command_unit, '(a)' ) 'set xlabel "T"' write ( command_unit, '(a)' ) 'set ylabel "D^2"' write ( command_unit, '(a)' ) & 'set title "Squared displacement (Red), Predicted (Black), Samples (Blue)"' write ( command_unit, '(a)' ) 'set grid' write ( command_unit, '(a)' ) 'set style data lines' write ( command_unit, '(a)' ) 'plot "' // trim ( data_filename ) // & '" using 1:2 title "sample 1" linecolor rgb "blue", \' write ( command_unit, '(a)' ) ' "' // trim ( data_filename ) // & '" using 1:3 title "sample 2" linecolor rgb "blue", \' write ( command_unit, '(a)' ) ' "' // trim ( data_filename ) // & '" using 1:4 title "sample 3" linecolor rgb "blue", \' write ( command_unit, '(a)' ) ' "' // trim ( data_filename ) // & '" using 1:5 title "sample 4" linecolor rgb "blue", \' write ( command_unit, '(a)' ) ' "' // trim ( data_filename ) // & '" using 1:6 title "sample 5" linecolor rgb "blue", \' write ( command_unit, '(a)' ) ' "' // trim ( data_filename ) // & '" using 1:7 title "Averaged" lw 3 linecolor rgb "red", \' write ( command_unit, '(a)' ) ' "' // trim ( data_filename ) // & '" using 1:8 title "Ideal" lw 3 linecolor rgb "black"' write ( command_unit, '(a)' ) 'quit' close ( unit = command_unit ) write ( *, '(a)' ) ' BROWNIAN_DISPLACEMENT plot commands stored in "' & // trim ( command_filename ) // '".' return end subroutine brownian_displacement_simulation ( k, n, m, d, t, seed, dsq ) !*****************************************************************************80 ! !! BROWNIAN_DISPLACEMENT_SIMULATION simulates Brownian displacement. ! ! Discussion: ! ! This function computes the square of the distance of the Brownian ! particle from the starting point, repeating this calculation ! several times. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 30 September 2012 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) K, the number of repetitions. ! ! Input, integer ( kind = 4 ) N, the number of time steps to take, plus 1. ! ! Input, integer ( kind = 4 ) M, the spatial dimension. ! ! Input, real ( kind = 8 ) D, the diffusion coefficient. ! Computationally, this is simply a scale factor between time and space. ! ! Input, real ( kind = 8 ) T, the total time. ! ! Input/output, integer ( kind = 4 ) SEED, a seed for the random ! number generator. ! ! Output, real ( kind = 8 ) DSQ(K,N), the displacements over time for each ! repetition. DSQ(:,1) is 0.0, because we include the displacement at the ! initial time. ! implicit none integer ( kind = 4 ) k integer ( kind = 4 ) m integer ( kind = 4 ) n real ( kind = 8 ) d real ( kind = 8 ) dsq(k,n) integer ( kind = 4 ) i integer ( kind = 4 ) seed real ( kind = 8 ) t real ( kind = 8 ) x(m,n) do i = 1, k call brownian_motion_simulation ( m, n, d, t, seed, x ) dsq(i,1:n) = sum ( x(1:m,1:n) ** 2, dim = 1 ) end do return end subroutine brownian_motion_display ( m, n, x, header ) !*****************************************************************************80 ! !! BROWNIAN_MOTION_DISPLAY displays successive Brownian motion positions. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 30 September 2012 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) M, the spatial dimension. ! M should be 1 or 2. ! ! Input, integer ( kind = 4 ) N, the number of time steps. ! ! Input, real ( kind = 8 ) X(M,N), the particle positions. ! ! Input, character ( len = * ) HEADER, an identifier for the output files. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n character ( len = 255 ) command_filename integer ( kind = 4 ) command_unit character ( len = 255 ) data_filename integer ( kind = 4 ) data_unit character ( len = * ) header integer ( kind = 4 ) i real ( kind = 8 ) t real ( kind = 8 ) x(m,n) if ( m /= 1 .and. m /= 2 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'BROWNIAN_MOTION_DISPLAY - Fatal error!' write ( *, '(a)' ) ' This routine can only handle M = 1 or 2.' stop end if ! ! Create the data file. ! call get_unit ( data_unit ) data_filename = trim ( header ) // '_data.txt' open ( unit = data_unit, file = data_filename, status = 'replace' ) if ( m == 1 ) then do i = 1, n t = real ( i - 1, kind = 8 ) / real ( n - 1, kind = 8 ) write ( data_unit, '(2x,g14.6,2x,g14.6)' ) t, x(1,i) end do else if ( m == 2 ) then do i = 1, n write ( data_unit, '(2x,g14.6,2x,g14.6)' ) x(1,i), x(2,i) end do end if close ( unit = data_unit ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' BROWNIAN_MOTION data stored in "' & // trim ( data_filename ) // '".' ! ! Create the command file. ! call get_unit ( command_unit ) command_filename = trim ( header ) // '_commands.txt' open ( unit = command_unit, file = command_filename, status = 'replace' ) write ( command_unit, '(a)' ) '# ' // trim ( command_filename ) write ( command_unit, '(a)' ) '#' write ( command_unit, '(a)' ) '# Usage:' write ( command_unit, '(a)' ) '# gnuplot < ' // trim ( command_filename ) write ( command_unit, '(a)' ) '#' write ( command_unit, '(a)' ) 'set term png' write ( command_unit, '(a)' ) 'set output "' // trim ( header ) // '.png"' write ( command_unit, '(a)' ) 'set xlabel "X"' write ( command_unit, '(a)' ) 'set ylabel "T"' write ( command_unit, '(a)' ) 'set title "Brownian motion in 1D"' write ( command_unit, '(a)' ) 'set grid' write ( command_unit, '(a)' ) 'set style data lines' write ( command_unit, '(a)' ) 'plot "' // trim ( data_filename ) & // '" using 1:2' write ( command_unit, '(a)' ) 'quit' close ( unit = command_unit ) write ( *, '(a)' ) ' BROWNIAN_MOTION plot commands stored in "' & // trim ( command_filename ) // '".' return end subroutine brownian_motion_simulation ( m, n, d, t, seed, x ) !*****************************************************************************80 ! !! BROWNIAN_MOTION_SIMULATION simulates Brownian motion. ! ! Discussion: ! ! Thanks to Feifei Xu for pointing out a missing factor of 2 in the ! stepsize calculation, 08 March 2016. ! ! Thanks to Joerg Peter Pfannmoeller for pointing out a missing factor ! of M in the stepsize calculation, 23 April 2018. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 23 April 2018 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) M, the spatial dimension. ! ! Input, integer ( kind = 4 ) N, the number of time steps to take, plus 1. ! ! Input, real ( kind = 8 ) D, the diffusion coefficient. ! ! Input, real ( kind = 8 ) T, the total time. ! ! Input/output, integer ( kind = 4 ) SEED, a seed for the random ! number generator. ! ! Output, real ( kind = 8 ) X(M,N), the initial position at time 0.0, and ! the N-1 successive locations of the particle. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n real ( kind = 8 ) d real ( kind = 8 ) dt real ( kind = 8 ) dx(m) integer ( kind = 4 ) j real ( kind = 8 ) norm_dx real ( kind = 8 ) r8_normal_01 real ( kind = 8 ) s integer ( kind = 4 ) seed real ( kind = 8 ) t real ( kind = 8 ) x(m,n) ! ! Set the time step. ! dt = t / real ( n - 1, kind = 8 ) ! ! Start at the origin. ! x(1:m,1) = 0.0D+00 ! ! Take N - 1 steps. ! do j = 2, n ! ! S is the stepsize. ! s = sqrt ( 2.0D+00 * real ( m, kind = 8 ) * d * dt ) * r8_normal_01 ( seed ) ! ! Direction DX is random, unit norm. ! if ( m == 1 ) then dx(1) = s else call r8vec_normal_01 ( m, seed, dx ) norm_dx = sqrt ( sum ( dx(1:m) ** 2 ) ) dx(1:m) = s * dx(1:m) / norm_dx end if ! ! Add the step to the current position. ! x(1:m,j) = x(1:m,j-1) + dx(1:m) 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 function i4_uniform ( a, b, seed ) !*****************************************************************************80 ! !! I4_UNIFORM returns a scaled pseudorandom I4. ! ! Discussion: ! ! An I4 is an integer ( kind = 4 ) value. ! ! The pseudorandom number will be scaled to be uniformly distributed ! between A and B. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 12 November 2006 ! ! Author: ! ! John Burkardt ! ! Reference: ! ! Paul Bratley, Bennett Fox, Linus Schrage, ! A Guide to Simulation, ! Springer Verlag, pages 201-202, 1983. ! ! Pierre L'Ecuyer, ! Random Number Generation, ! in Handbook of Simulation, ! edited by Jerry Banks, ! Wiley Interscience, page 95, 1998. ! ! 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 ) A, B, the limits of the interval. ! ! Input/output, integer ( kind = 4 ) SEED, the "seed" value, which ! should NOT be 0. On output, SEED has been updated. ! ! Output, integer ( kind = 4 ) I4_UNIFORM, a number between A and B. ! implicit none integer ( kind = 4 ) a integer ( kind = 4 ) b integer ( kind = 4 ), parameter :: i4_huge = 2147483647 integer ( kind = 4 ) i4_uniform integer ( kind = 4 ) k real ( kind = 4 ) r integer ( kind = 4 ) seed integer ( kind = 4 ) value if ( seed == 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'I4_UNIFORM - Fatal error!' write ( *, '(a)' ) ' Input value of SEED = 0.' stop end if k = seed / 127773 seed = 16807 * ( seed - k * 127773 ) - k * 2836 if ( seed < 0 ) then seed = seed + i4_huge end if r = real ( seed, kind = 4 ) * 4.656612875E-10 ! ! Scale R to lie between A-0.5 and B+0.5. ! r = ( 1.0E+00 - r ) * ( real ( min ( a, b ), kind = 4 ) - 0.5E+00 ) & + r * ( real ( max ( a, b ), kind = 4 ) + 0.5E+00 ) ! ! Use rounding to convert R to an integer between A and B. ! value = nint ( r, kind = 4 ) value = max ( value, min ( a, b ) ) value = min ( value, max ( a, b ) ) i4_uniform = value return end function r8_normal_01 ( seed ) !*****************************************************************************80 ! !! R8_NORMAL_01 returns a unit pseudonormal R8. ! ! Discussion: ! ! The standard normal probability distribution function (PDF) has ! mean 0 and standard deviation 1. ! ! Because this routine uses the Box Muller method, it requires pairs ! of uniform random values to generate a pair of normal random values. ! This means that on every other call, essentially, the input value of ! SEED is ignored, since the code saves the second normal random value. ! ! If you didn't know this, you might be confused since, usually, the ! output of a random number generator can be completely controlled by ! the input value of the SEED. If I were more careful, I could rewrite ! this routine so that it would distinguish between cases where the input ! value of SEED is the output value from the previous call (all is well) ! and those cases where it is not (the user has decided to do something ! new. Restart the uniform random number sequence.) But I'll leave ! that for later. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 17 July 2006 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input/output, integer ( kind = 4 ) SEED, a seed for the random number ! generator. ! ! Output, real ( kind = 8 ) R8_NORMAL_01, a sample of the standard ! normal PDF. ! implicit none real ( kind = 8 ), parameter :: pi = 3.141592653589793D+00 real ( kind = 8 ) r1 real ( kind = 8 ) r2 real ( kind = 8 ) r8_normal_01 real ( kind = 8 ) r8_uniform_01 integer ( kind = 4 ) seed integer ( kind = 4 ), save :: seed2 = 0 integer ( kind = 4 ), save :: used = 0 real ( kind = 8 ) x real ( kind = 8 ), save :: y = 0.0D+00 ! ! On odd numbered calls, generate two uniforms, create two normals, ! return the first normal and its corresponding seed. ! if ( mod ( used, 2 ) == 0 ) then r1 = r8_uniform_01 ( seed ) if ( r1 == 0.0D+00 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8_NORMAL_01 - Fatal error!' write ( *, '(a)' ) ' R8_UNIFORM_01 returned a value of 0.' stop end if seed2 = seed r2 = r8_uniform_01 ( seed2 ) x = sqrt ( - 2.0D+00 * log ( r1 ) ) * cos ( 2.0D+00 * pi * r2 ) y = sqrt ( - 2.0D+00 * log ( r1 ) ) * sin ( 2.0D+00 * pi * r2 ) ! ! On odd calls, return the second normal and its corresponding seed. ! else seed = seed2 x = y end if used = used + 1 r8_normal_01 = x return end function r8_uniform_01 ( seed ) !*****************************************************************************80 ! !! R8_UNIFORM_01 returns a unit pseudorandom R8. ! ! Discussion: ! ! An R8 is a real ( kind = 8 ) value. ! ! For now, the input quantity SEED is an integer variable. ! ! This routine implements the recursion ! ! seed = 16807 * seed mod ( 2^31 - 1 ) ! r8_uniform_01 = seed / ( 2^31 - 1 ) ! ! The integer arithmetic never requires more than 32 bits, ! including a sign bit. ! ! If the initial seed is 12345, then the first three computations are ! ! Input Output R8_UNIFORM_01 ! SEED SEED ! ! 12345 207482415 0.096616 ! 207482415 1790989824 0.833995 ! 1790989824 2035175616 0.947702 ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 05 July 2006 ! ! Author: ! ! John Burkardt ! ! Reference: ! ! Paul Bratley, Bennett Fox, Linus Schrage, ! A Guide to Simulation, ! Springer Verlag, pages 201-202, 1983. ! ! Pierre L'Ecuyer, ! Random Number Generation, ! in Handbook of Simulation, ! edited by Jerry Banks, ! Wiley Interscience, page 95, 1998. ! ! 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/output, integer ( kind = 4 ) SEED, the "seed" value, which should ! NOT be 0. On output, SEED has been updated. ! ! Output, real ( kind = 8 ) R8_UNIFORM_01, a new pseudorandom variate, ! strictly between 0 and 1. ! implicit none integer ( kind = 4 ), parameter :: i4_huge = 2147483647 integer ( kind = 4 ) k real ( kind = 8 ) r8_uniform_01 integer ( kind = 4 ) seed if ( seed == 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8_UNIFORM_01 - Fatal error!' write ( *, '(a)' ) ' Input value of SEED = 0.' stop end if k = seed / 127773 seed = 16807 * ( seed - k * 127773 ) - k * 2836 if ( seed < 0 ) then seed = seed + i4_huge end if r8_uniform_01 = real ( seed, kind = 8 ) * 4.656612875D-10 return end subroutine r8vec_normal_01 ( n, seed, x ) !*****************************************************************************80 ! !! R8VEC_NORMAL_01 returns a unit pseudonormal R8VEC. ! ! Discussion: ! ! An R8VEC is a vector of R8's. ! ! The standard normal probability distribution function (PDF) has ! mean 0 and standard deviation 1. ! ! This routine can generate a vector of values on one call. It ! has the feature that it should provide the same results ! in the same order no matter how we break up the task. ! ! Before calling this routine, the user may call RANDOM_SEED ! in order to set the seed of the random number generator. ! ! The Box-Muller method is used, which is efficient, but ! generates an even number of values each time. On any call ! to this routine, an even number of new values are generated. ! Depending on the situation, one value may be left over. ! In that case, it is saved for the next call. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 17 July 2006 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, the number of values desired. If N is ! negative,then the code will flush its internal memory; in particular, ! if there is a saved value to be used on the next call, it is ! instead discarded. This is useful if the user has reset the ! random number seed, for instance. ! ! Input/output, integer ( kind = 4 ) SEED, a seed for the random ! number generator. ! ! Output, real ( kind = 8 ) X(N), a sample of the standard normal PDF. ! ! Local parameters: ! ! Local, integer ( kind = 4 ) MADE, records the number of values that have ! been computed. On input with negative N, this value overwrites ! the return value of N, so the user can get an accounting of ! how much work has been done. ! ! Local, real ( kind = 8 ) R(N+1), is used to store some uniform ! random values. Its dimension is N+1, but really it is only needed ! to be the smallest even number greater than or equal to N. ! ! Local, integer SAVED, is 0 or 1 depending on whether there is a ! single saved value left over from the previous call. ! ! Local, integer X_LO_INDEX, X_HI_INDEX, records the range of entries of ! X that we need to compute. This starts off as 1:N, but is adjusted ! if we have a saved value that can be immediately stored in X(1), ! and so on. ! ! Local, real ( kind = 8 ) Y, the value saved from the previous call, if ! SAVED is 1. ! implicit none integer ( kind = 4 ) n integer ( kind = 4 ) m integer ( kind = 4 ), save :: made = 0 real ( kind = 8 ), parameter :: pi = 3.141592653589793D+00 real ( kind = 8 ) r(n+1) real ( kind = 8 ) r8_uniform_01 integer ( kind = 4 ), save :: saved = 0 integer ( kind = 4 ) seed real ( kind = 8 ) x(n) integer ( kind = 4 ) x_hi_index integer ( kind = 4 ) x_lo_index real ( kind = 8 ), save :: y = 0.0D+00 ! ! I'd like to allow the user to reset the internal data. ! But this won't work properly if we have a saved value Y. ! I'm making a crock option that allows the user to signal ! explicitly that any internal memory should be flushed, ! by passing in a negative value for N. ! if ( n < 0 ) then n = made made = 0 saved = 0 y = 0.0D+00 return else if ( n == 0 ) then return end if ! ! Record the range of X we need to fill in. ! x_lo_index = 1 x_hi_index = n ! ! Use up the old value, if we have it. ! if ( saved == 1 ) then x(1) = y saved = 0 x_lo_index = 2 end if ! ! Maybe we don't need any more values. ! if ( x_hi_index - x_lo_index + 1 == 0 ) then ! ! If we need just one new value, do that here to avoid null arrays. ! else if ( x_hi_index - x_lo_index + 1 == 1 ) then r(1) = r8_uniform_01 ( seed ) if ( r(1) == 0.0D+00 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8VEC_NORMAL_01 - Fatal error!' write ( *, '(a)' ) ' R8_UNIFORM_01 returned a value of 0.' stop end if r(2) = r8_uniform_01 ( seed ) x(x_hi_index) = & sqrt ( -2.0D+00 * log ( r(1) ) ) * cos ( 2.0D+00 * pi * r(2) ) y = sqrt ( -2.0D+00 * log ( r(1) ) ) * sin ( 2.0D+00 * pi * r(2) ) saved = 1 made = made + 2 ! ! If we require an even number of values, that's easy. ! else if ( mod ( x_hi_index - x_lo_index + 1, 2 ) == 0 ) then m = ( x_hi_index - x_lo_index + 1 ) / 2 call r8vec_uniform_01 ( 2*m, seed, r ) x(x_lo_index:x_hi_index-1:2) = & sqrt ( -2.0D+00 * log ( r(1:2*m-1:2) ) ) & * cos ( 2.0D+00 * pi * r(2:2*m:2) ) x(x_lo_index+1:x_hi_index:2) = & sqrt ( -2.0D+00 * log ( r(1:2*m-1:2) ) ) & * sin ( 2.0D+00 * pi * r(2:2*m:2) ) made = made + x_hi_index - x_lo_index + 1 ! ! If we require an odd number of values, we generate an even number, ! and handle the last pair specially, storing one in X(N), and ! saving the other for later. ! else x_hi_index = x_hi_index - 1 m = ( x_hi_index - x_lo_index + 1 ) / 2 + 1 call r8vec_uniform_01 ( 2*m, seed, r ) x(x_lo_index:x_hi_index-1:2) = & sqrt ( -2.0D+00 * log ( r(1:2*m-3:2) ) ) & * cos ( 2.0D+00 * pi * r(2:2*m-2:2) ) x(x_lo_index+1:x_hi_index:2) = & sqrt ( -2.0D+00 * log ( r(1:2*m-3:2) ) ) & * sin ( 2.0D+00 * pi * r(2:2*m-2:2) ) x(n) = sqrt ( -2.0D+00 * log ( r(2*m-1) ) ) & * cos ( 2.0D+00 * pi * r(2*m) ) y = sqrt ( -2.0D+00 * log ( r(2*m-1) ) ) & * sin ( 2.0D+00 * pi * r(2*m) ) saved = 1 made = made + x_hi_index - x_lo_index + 2 end if return end subroutine r8vec_uniform_01 ( n, seed, r ) !*****************************************************************************80 ! !! R8VEC_UNIFORM_01 returns a unit pseudorandom R8VEC. ! ! Discussion: ! ! An R8VEC is a vector of R8's. ! ! For now, the input quantity SEED is an integer ( kind = 4 ) variable. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 05 July 2006 ! ! 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 ) N, the number of entries in the vector. ! ! 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(N), the vector of pseudorandom values. ! implicit none integer ( kind = 4 ) n integer ( kind = 4 ) i integer ( kind = 4 ) k integer ( kind = 4 ) seed real ( kind = 8 ) r(n) if ( seed == 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8VEC_UNIFORM_01 - Fatal error!' write ( *, '(a)' ) ' Input value of SEED = 0.' stop end if do i = 1, n k = seed / 127773 seed = 16807 * ( seed - k * 127773 ) - k * 2836 if ( seed < 0 ) then seed = seed + 2147483647 end if r(i) = real ( seed, kind = 8 ) * 4.656612875D-10 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