program main c*****************************************************************************80 c cc trapezoidal_test() tests trapezoidal(). c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 30 November 2023 c c Author: c c John Burkardt c implicit none integer n double precision tspan(2) double precision y0(2) call timestamp ( ) write ( *, '(a)' ) '' write ( *, '(a)' ) 'trapezoidal_test():' write ( *, '(a)' ) ' Fortran77 version' write ( *, '(a)' ) ' Test trapezoidal().' tspan(1) = 0.0D+00 tspan(2) = 10.0D+00 y0(1) = 5000.0D+00 y0(2) = 100.0D+00 n = 1000 call predator_prey_tr_test ( tspan, y0, n ) tspan(1) = 0.0D+00 tspan(2) = 1.0D+00 y0(1) = 0.0D+00 n = 27 call stiff_tr_test ( tspan, y0(1), n ) c c Terminate. c write ( *, '(a)' ) '' write ( *, '(a)' ) 'trapezoidal_test():' write ( *, '(a)' ) ' Normal end of execution.' write ( *, '(a)' ) '' call timestamp ( ) stop 0 end subroutine predator_prey_tr_test ( tspan, y0, n ) c*********************************************************************72 c cc predator_prey_tr_test(): predator-prey using trapezoidal(). c c Discussion: c c The physical system under consideration is a pair of animal populations. c c The PREY reproduce rapidly for each animal alive at the beginning of the c year, two more will be born by the end of the year. The prey do not have c a natural death rate instead, they only die by being eaten by the predator. c Every prey animal has 1 chance in 1000 of being eaten in a given year by c a given predator. c c The PREDATORS only die of starvation, but this happens very quickly. c If unfed, a predator will tend to starve in about 1/10 of a year. c On the other hand, the predator reproduction rate is dependent on c eating prey, and the chances of this depend on the number of available prey. c c The resulting differential equations can be written: c c PREY(0) = 5000 c PRED(0) = 100 c c d PREY / dT = 2 * PREY(T) - 0.001 * PREY(T) * PRED(T) c d PRED / dT = - 10 * PRED(T) + 0.002 * PREY(T) * PRED(T) c c Here, the initial values (5000,100) are a somewhat arbitrary starting point. c c The pair of ordinary differential equations that result have an interesting c behavior. For certain choices of the interaction coefficients (such as c those given here), the populations of predator and prey will tend to c a periodic oscillation. The two populations will be out of phase the number c of prey will rise, then after a delay, the predators will rise as the prey c begins to fall, causing the predator population to crash again. c c There is a conserved quantity, which here would be: c E(r,f) = 0.002 r + 0.001 f - 10 ln(r) - 2 ln(f) c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 30 November 2023 c c Author: c c John Burkardt c c Reference: c c George Lindfield, John Penny, c Numerical Methods Using MATLAB, c Second Edition, c Prentice Hall, 1999, c ISBN: 0-13-012641-1, c LC: QA297.P45. c c Input: c c double precision tspan: contains [ T0, TMAX ], the initial and final c times. A reasonable value might be [ 0, 5 ]. c c double precision y0 = [ PREY, PRED ], the initial number of prey and c predators. A reasonable value might be [ 5000, 100 ]. c c integer n: the number of time steps to take. c implicit none integer, parameter :: m = 2 integer n character ( len = 255 ) command_filename integer command_unit character ( len = 255 ) data_filename integer data_unit character ( len = * ), parameter :: header = 'predator_prey' integer i external predator_prey_dydt double precision t(n+1) double precision tspan(2) double precision y0(m) double precision y(n+1,m) write ( *, '(a)' ) '' write ( *, '(a)' ) 'predator_prey_tr_test()' write ( *, '(a)' ) '' write ( *, '(a)' ) ' A pair of ordinary differential equations' write ( *, '(a)' ) ' for a population of predators and prey ' write ( *, '(a)' ) ' are solved using trapezoidal().' write ( *, '(a)' ) '' write ( *, '(a)' ) ' The exact solution shows periodic behavior' write ( *, '(a)' ) ' with a fixed period and amplitude.' call trapezoidal ( predator_prey_dydt, tspan, y0, n, m, t, y ) c c Create the data file. c call get_unit ( data_unit ) data_filename = header // '_data.txt' open ( unit = data_unit, file = data_filename, & status = 'replace' ) do i = 1, n + 1 write ( data_unit, '(5(2x,g14.6))' ) t(i), y(i,1), y(i,2) end do close ( unit = data_unit ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' predator_prey_tr_test: data stored in "' & // trim ( data_filename ) // '".' c c Create the command file. c 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 "<-- Prey -->"' write ( command_unit, '(a)' ) 'set ylabel "<-- Predator -->"' write ( command_unit, '(a)' ) & 'set title "predator prey ODE, trapezoidal"' write ( command_unit, '(a)' ) 'set grid' write ( command_unit, '(a)' ) 'set style data lines' write ( command_unit, '(a)' ) & 'plot "' // trim ( data_filename ) // & '" using 2:3 with lines lw 3' write ( command_unit, '(a)' ) 'quit' close ( unit = command_unit ) write ( *, '(a)' ) & ' predator_prey_tr_test: plot commands stored in "' & // trim ( command_filename ) // '".' return end subroutine predator_prey_dydt ( t, y, dydt ) c*********************************************************************72 c cc predator_prey_dydt() evaluates the right hand side of the system. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 10 November 2023 c c Author: c c John Burkardt c c Reference: c c George Lindfield, John Penny, c Numerical Methods Using MATLAB, c Second Edition, c Prentice Hall, 1999, c ISBN: 0-13-012641-1, c LC: QA297.P45. c c Input: c c double precision t: the current time. c c double precision y(2): the current solution variables, rabbits and foxes. c c Output: c c double precision dydt(2): the right hand side of the 2 ODE's. c implicit none double precision dydt(2) double precision t double precision y(2) call r8_fake_use ( t ) dydt(1) = 2.0 * y(1) - 0.001 * y(1) * y(2) dydt(2) = - 10.0 * y(2) + 0.002 * y(1) * y(2) return end subroutine stiff_tr_test ( tspan, y0, n ) c*********************************************************************72 c cc stiff_tr_test() uses the trapezoidal method on the stiff ODE. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 30 November 2023 c c Author: c c John Burkardt c c Input: c c double precision tspan(2): the first and last times. c c double precision Y0: the initial condition. c c integer N: the number of steps to take. c implicit none integer, parameter :: m = 1 integer n integer, parameter :: n2 = 101 external stiff_dydt double precision t1(n+1) double precision t2(n2) double precision tspan(2) double precision y1(n+1) double precision y0 double precision y2(n2) write ( *, '(a)' ) '' write ( *, '(a)' ) 'stiff_tr_test():' write ( *, '(a)' ) ' Solve stiff ODE using backward_euler()' call trapezoidal ( stiff_dydt, tspan, y0, n, m, t1, y1 ) call r8vec_linspace ( n2, tspan(1), tspan(2), t2 ) call stiff_exact ( n2, t2, y2 ) call plot2 ( n+1, t1, y1, n2, t2, y2, 'stiff', & 'stiff (trapezoidal)' ) return end subroutine stiff_dydt ( t, y, dydt ) c*********************************************************************72 c cc stiff_dydt() evaluates the right hand side of the stiff ODE. c c Discussion: c c y' = 50 * ( cos(t) - y ) c y(0) = 0 c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 10 November 2023 c c Author: c c John Burkardt c c Input: c c double precision T, Y(1): the time and solution value. c c Output: c c double precision DYDT(1): the derivative value. c implicit none double precision dydt(1) double precision t double precision y(1) dydt(1) = 50.0D+00 * ( cos ( t ) - y(1) ) return end subroutine stiff_exact ( n, t, y ) c*********************************************************************72 c cc stiff_exact() evaluates the exact solution of the stiff ODE. c c Discussion: c c y' = 50 * ( cos(t) - y ) c y(0) = 0 c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 27 April 2020 c c Author: c c John Burkardt c c Input: c c integer N: the number of values. c c double precision T(N): the evaluation times. c c Output: c c double precision Y(N): the exact solution values. c implicit none integer n integer i double precision t(n) double precision y(n) do i = 1, n y(i) = 50.0D+00 * ( sin ( t(i) ) + 50.0D+00 * cos ( t(i) ) & - 50.0D+00 * exp ( - 50.0D+00 * t(i) ) ) / 2501.0D+00 end do return end subroutine get_unit ( iunit ) c*********************************************************************72 c cc get_unit() returns a free FORTRAN unit number. c c Discussion: c c A "free" FORTRAN unit number is a value between 1 and 99 which c is not currently associated with an I/O device. A free FORTRAN unit c number is needed in order to open a file with the OPEN command. c c If IUNIT = 0, then no free FORTRAN unit could be found, although c all 99 units were checked (except for units 5, 6 and 9, which c are commonly reserved for console I/O). c c Otherwise, IUNIT is a value between 1 and 99, representing a c free FORTRAN unit. Note that GET_UNIT assumes that units 5 and 6 c are special, and will never return those values. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 02 September 2013 c c Author: c c John Burkardt c c Parameters: c c Output, integer IUNIT, the free unit number. c implicit none integer i integer iunit logical value iunit = 0 do i = 1, 99 if ( i .ne. 5 .and. i .ne. 6 .and. i .ne. 9 ) then inquire ( unit = i, opened = value, err = 10 ) if ( .not. value ) then iunit = i return end if end if 10 continue end do return end subroutine plot2 ( n1, t1, y1, n2, t2, y2, header, title ) c*********************************************************************72 c cc plot2() plots two curves together. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 27 April 2020 c c Author: c c John Burkardt c c Input: c c integer N1: the size of the first data set. c c double precision T1(N1), Y1(N1), the first dataset. c c integer N2: the size of the second data set. c c double precision T2(N2), Y2(N2), the secod dataset. c c character ( len = * ) HEADER: an identifier for the data. c c character ( len = * ) TITLE: a title to appear in the plot. c implicit none integer, parameter :: rk = kind ( 1.0D+00 ) integer n1 integer n2 character ( len = 255 ) command_filename integer command_unit character ( len = 255 ) data1_filename character ( len = 255 ) data2_filename integer data_unit character ( len = * ) header integer i character ( len = * ) title double precision t1(n1) double precision t2(n2) double precision y1(n1) double precision y2(n2) c c Create the data files. c call get_unit ( data_unit ) data1_filename = header // '_data1.txt' open ( unit = data_unit, file = data1_filename, & status = 'replace' ) do i = 1, n1 write ( data_unit, '(5(2x,g14.6))' ) t1(i), y1(i) end do close ( unit = data_unit ) call get_unit ( data_unit ) data2_filename = header // '_data2.txt' open ( unit = data_unit, file = data2_filename, & status = 'replace' ) do i = 1, n2 write ( data_unit, '(5(2x,g14.6))' ) t2(i), y2(i) end do close ( unit = data_unit ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' plot2: data stored in "' & // trim ( data1_filename ) // '" and "' & // trim ( data2_filename ) // '".' c c Create the command file. c 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 "<-- Y(T) -->"' write ( command_unit, '(a)' ) 'set title "' & // trim ( title ) // '"' write ( command_unit, '(a)' ) 'set grid' write ( command_unit, '(a)' ) 'set style data lines' write ( command_unit, '(a)' ) 'plot "' & // trim ( data1_filename ) // & '" using 1:2 with lines lw 3 lt rgb "red",\' write ( command_unit, '(a)' ) ' "' & // trim ( data2_filename ) // & '" using 1:2 with lines lw 3 lt rgb "blue"' write ( command_unit, '(a)' ) 'quit' close ( unit = command_unit ) write ( *, '(a)' ) ' plot2: plot commands stored in "' & // trim ( command_filename ) // '".' return end subroutine r8_fake_use ( x ) c*********************************************************************72 c cc r8_fake_use() pretends to use an R8 variable. c c Discussion: c c Some compilers will issue a warning if a variable is unused. c Sometimes there's a good reason to include a variable in a program, c but not to use it. Calling this function with that variable as c the argument will shut the compiler up. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 19 August 2023 c c Author: c c John Burkardt c c Input: c c double precision x: the variable to be "used". c implicit none double precision x if ( x /= x ) then write ( *, '(a)' ) ' r8_fake_use(): variable is NAN.' end if return end subroutine r8vec_linspace ( n, a, b, x ) c*********************************************************************72 c cc r8vec_linspace() creates a vector of linearly spaced values. c c Discussion: c c An R8VEC is a vector of R8's. c c 4 points evenly spaced between 0 and 12 will yield 0, 4, 8, 12. c c In other words, the interval is divided into N-1 even subintervals, c and the endpoints of intervals are used as the points. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 14 March 2011 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the number of entries in the vector. c c Input, double precision A, B, the first and last entries. c c Output, double precision X(N), a vector of linearly spaced data. c implicit none integer n double precision a double precision b integer i double precision x(n) if ( n .eq. 1 ) then x(1) = ( a + b ) / 2.0D+00 else do i = 1, n x(i) = ( dble ( n - i ) * a & + dble ( i - 1 ) * b ) & / dble ( n - 1 ) end do end if return end subroutine timestamp ( ) c*********************************************************************72 c cc timestamp() prints the YMDHMS date as a timestamp. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 12 June 2014 c c Author: c c John Burkardt c implicit none character * ( 8 ) ampm integer d character * ( 8 ) date integer h integer m integer mm character * ( 9 ) month(12) integer n integer s character * ( 10 ) time integer y save month data month / & 'January ', 'February ', 'March ', 'April ', & 'May ', 'June ', 'July ', 'August ', & 'September', 'October ', 'November ', 'December ' / call date_and_time ( date, time ) read ( date, '(i4,i2,i2)' ) y, m, d read ( time, '(i2,i2,i2,1x,i3)' ) h, n, s, mm if ( h .lt. 12 ) then ampm = 'AM' else if ( h .eq. 12 ) then if ( n .eq. 0 .and. s .eq. 0 ) then ampm = 'Noon' else ampm = 'PM' end if else h = h - 12 if ( h .lt. 12 ) then ampm = 'PM' else if ( h .eq. 12 ) then if ( n .eq. 0 .and. s .eq. 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