program main !*****************************************************************************80 ! !! heated_plate() solves the steady heat equation on a rectangle. ! ! Discussion: ! ! This code solves the steady state heat equation on a rectangular region. ! ! The sequential version of this program needs approximately ! 18/eps iterations to complete. ! ! ! The physical region, and the boundary conditions, are suggested ! by this diagram; ! ! W = 0 ! +------------------+ ! | | ! W = 100 | | W = 100 ! | | ! +------------------+ ! W = 100 ! ! The region is covered with a grid of M by N nodes, and an N by N ! array W is used to record the temperature. The correspondence between ! array indices and locations in the region is suggested by giving the ! indices of the four corners: ! ! I = 0 ! [0][0]-------------[0][N-1] ! | | ! J = 0 | | J = N-1 ! | | ! [M-1][0]-----------[M-1][N-1] ! I = M-1 ! ! The steady state solution to the discrete heat equation satisfies the ! following condition at an interior grid point: ! ! W[Central] = (1/4) * ( W[North] + W[South] + W[East] + W[West] ) ! ! where "Central" is the index of the grid point, "North" is the index ! of its immediate neighbor to the "north", and so on. ! ! Given an approximate solution of the steady state heat equation, a ! "better" solution is given by replacing each interior point by the ! average of its 4 neighbors - in other words, by using the condition ! as an ASSIGNMENT statement: ! ! W[Central] <= (1/4) * ( W[North] + W[South] + W[East] + W[West] ) ! ! If this process is repeated often enough, the difference between ! successive estimates of the solution will go to zero. ! ! This program carries out such an iteration, using a tolerance specified by ! the user, and writes the final estimate of the solution to a file that can ! be used for graphic processing. ! ! Thanks to Vivek Rao for updating the command line argument handling, ! 20 February 2025. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 02 April 2025 ! ! Author: ! ! Original FORTRAN version by Michael Quinn. ! This version by John Burkardt. ! ! Reference: ! ! Michael Quinn, ! Parallel Programming in C with MPI and OpenMP, ! McGraw-Hill, 2004, ! ISBN13: 978-0071232654, ! LC: QA76.73.C15.Q55. ! ! Parameters: ! ! real ( kind = rk8 ) EPS, the error tolerance. ! ! Local: ! ! real ( kind = rk8 ) DIFF, the norm of the change in the solution from ! one iteration to the next. ! ! real ( kind = rk8 ) MEAN, the average of the boundary values, used ! to initialize the values of the solution in the interior. ! ! real ( kind = rk8 ) U(M,N), the solution at the previous iteration. ! ! real ( kind = rk8 ) W(M,N), the solution computed at the latest ! iteration. ! implicit none integer, parameter :: rk8 = kind ( 1.0D+00 ) integer, parameter :: m = 200 integer, parameter :: n = 200 real ( kind = rk8 ) ctime real ( kind = rk8 ) ctime1 real ( kind = rk8 ) ctime2 real ( kind = rk8 ) diff real ( kind = rk8 ) eps character ( len = 255 ) header integer iterations integer iterations_print real ( kind = rk8 ) mean real ( kind = rk8 ), allocatable :: u(:,:) real ( kind = rk8 ), allocatable :: w(:,:) real ( kind = rk8 ), allocatable :: x(:) real ( kind = rk8 ), allocatable :: y(:) call timestamp ( ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'heated_plate():' write ( *, '(a)' ) ' Fortran90 version' write ( *, '(a)' ) & ' A program to solve for the steady state temperature distribution' write ( *, '(a)' ) ' over a rectangular plate.' write ( *, '(a)' ) ' ' write ( *, '(a,i8,a,i8,a)' ) ' Spatial grid of ', m, ' by ', n, ' points.' eps = 0.0001D+00 write ( *, '(a)' ) ' ' write ( *, '(a,g14.6)' ) & ' The iteration will repeat until the change is <= ', eps diff = eps ! ! Set X and Y. ! allocate ( x(1:n) ) allocate ( y(1:m) ) call r8vec_linspace ( n, 0.0D+00, 1.0D+00, x ) call r8vec_linspace ( m, 0.0D+00, 1.0D+00, y ) ! ! Set the boundary values, which don't change. ! allocate ( u(1:m,1:n) ) allocate ( w(1:m,1:n) ) w(2:m-1,1) = 100.0 w(2:m-1,n) = 100.0 w(m,1:n) = 100.0 w(1,1:n) = 0.0 ! ! Average the boundary values, to come up with a reasonable ! initial value for the interior. ! mean = ( & sum ( w(2:m-1,1) ) & + sum ( w(2:m-1,n) ) & + sum ( w(m,1:n) ) & + sum ( w(1,1:n) ) ) & / dble ( 2 * m + 2 * n - 4 ) ! ! Initialize the interior solution to the mean value. ! w(2:m-1,2:n-1) = mean ! ! iterate until the new solution W differs from the old solution U ! by no more than EPS. ! iterations = 0 iterations_print = 1 write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Iteration Change' write ( *, '(a)' ) ' ' call cpu_time ( ctime1 ) do while ( eps <= diff ) u(1:m,1:n) = w(1:m,1:n) w(2:m-1,2:n-1) = 0.25 * ( & u(1:m-2,2:n-1) & + u(3:m,2:n-1) & + u(2:m-1,1:n-2) & + u(2:m-1,3:n) ) diff = maxval ( abs ( u(1:m,1:n) - w(1:m,1:n) ) ) iterations = iterations + 1 if ( iterations == iterations_print ) then write ( *, '(2x,i8,2x,g14.6)' ) iterations, diff iterations_print = 2 * iterations_print end if end do call cpu_time ( ctime2 ) ctime = ctime2 - ctime write ( *, '(a)' ) ' ' write ( *, '(2x,i8,2x,g14.6)' ) iterations, diff write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Error tolerance achieved.' write ( *, '(a,g14.6)' ) ' CPU time = ', ctime ! ! Write the solution to a file. ! GNUPLOT requires a regular grid, and a blank line before a new ! Y value. ! header = 'heated_plate' call gnuplot_fxy ( m, n, x, y, w, trim ( header ) ) ! ! Terminate. ! write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'heated_plate():' write ( *, '(a)' ) ' Normal end of execution.' write ( *, '(a)' ) '' call timestamp ( ) deallocate ( u ) deallocate ( w ) deallocate ( x ) deallocate ( y ) stop 0 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. The code assumes that units 5 and 6 ! are special, and will never return those values. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 26 October 2008 ! ! Author: ! ! John Burkardt ! ! 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 gnuplot_fxy ( m, n, x, y, z, header ) !*****************************************************************************80 ! !! gnuplot_fxy() plots z(x,y) data. ! ! Discussion: ! ! Actually, we simply create two files for processing by gnuplot(). ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 02 April 2025 ! ! Author: ! ! John Burkardt ! implicit none integer, parameter :: rk8 = kind ( 1.0D+00 ) integer m integer n character ( len = 255 ) command_filename integer command_unit character ( len = 255 ) data_filename integer data_unit character ( len = * ) header integer i integer j character ( len = 255 ) png_filename real ( kind = rk8 ) x(n) real ( kind = rk8 ) y(m) real ( kind = rk8 ) z(m,n) ! ! Create a graphics data file. ! data_filename = header // '_data.txt' call get_unit ( data_unit ) open ( unit = data_unit, file = data_filename, status = 'replace' ) do i = 1, m do j = 1, n write ( data_unit, '(3g14.6)' ) x(j), y(i), z(i,j) end do if ( i < m ) then write ( data_unit, '(a)' ) '' end if end do close ( unit = data_unit ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Created graphics data file "' & // trim ( data_filename ) // '".' png_filename = header // '.png' ! ! Create graphics command file. ! call get_unit ( command_unit ) command_filename = 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 ( png_filename ) // '"' write ( command_unit, '(a)' ) 'set xlabel "<--- X --->"' write ( command_unit, '(a)' ) 'set ylabel "<--- Y --->"' write ( command_unit, '(a)' ) 'set title "' // header // '"' write ( command_unit, '(a)' ) 'set grid' write ( command_unit, '(a)' ) 'unset surface' write ( command_unit, '(a)' ) 'set contour base' write ( command_unit, '(a)' ) 'set view map' write ( command_unit, '(a)' ) 'set pm3d' write ( command_unit, '(a)' ) 'splot "' // trim ( data_filename ) // & '" with pm3d' close ( unit = command_unit ) write ( *, '(a)' ) & ' Created command file "' // trim ( command_filename ) // '".' return end subroutine r8vec_linspace ( n, a, b, x ) !*****************************************************************************80 ! !! r8vec_linspace() creates a vector of linearly spaced values. ! ! Discussion: ! ! An R8VEC is a vector of R8's. ! ! 4 points evenly spaced between 0 and 12 will yield 0, 4, 8, 12. ! ! In other words, the interval is divided into N-1 even subintervals, ! and the endpoints of intervals are used as the points. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 14 March 2011 ! ! Author: ! ! John Burkardt ! ! Input: ! ! integer N, the number of entries in the vector. ! ! real ( kind = rk8 ) A, B, the first and last entries. ! ! Output: ! ! real ( kind = rk8 ) X(N), a vector of linearly spaced data. ! implicit none integer, parameter :: rk8 = kind ( 1.0D+00 ) integer n real ( kind = rk8 ) a real ( kind = rk8 ) b integer i real ( kind = rk8 ) x(n) if ( n == 1 ) then x(1) = ( a + b ) / 2.0D+00 else do i = 1, n x(i) = ( real ( n - i, kind = rk8 ) * a & + real ( i - 1, kind = rk8 ) * b ) & / real ( n - 1, kind = rk8 ) end do end if 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.2,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