PartMC 2.1.2
util.F90
Go to the documentation of this file.
00001 ! Copyright (C) 2005-2011 Nicole Riemer and Matthew West
00002 ! Licensed under the GNU General Public License version 2 or (at your
00003 ! option) any later version. See the file COPYING for details.
00004 
00005 !> \file
00006 !> The pmc_util module.
00007 
00008 !> Common utility subroutines.
00009 module pmc_util
00010 
00011   use pmc_constants
00012 #ifdef PMC_USE_MPI
00013   use mpi
00014 #endif
00015 
00016   !> Maximum number of IO units usable simultaneously.
00017   integer, parameter :: max_units = 200
00018   !> Minimum unit number to allocate.
00019   integer, parameter :: unit_offset = 10
00020   !> Table of unit numbers storing allocation status.
00021   logical, save :: unit_used(max_units) = .false.
00022 
00023   !> Length of string for converting numbers.
00024   integer, parameter :: PMC_UTIL_CONVERT_STRING_LEN = 100
00025   !> Maximum length of filenames.
00026   integer, parameter :: PMC_MAX_FILENAME_LEN = 300
00027 
00028 contains
00029   
00030 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00031 
00032   !> Prints a warning message.
00033   subroutine warn_msg(code, warning_msg)
00034 
00035     !> Status code to use.
00036     integer, intent(in) :: code
00037     !> Message to display.
00038     character(len=*), intent(in) :: warning_msg
00039 
00040     write(0,'(a)') 'WARNING (PartMC-' // trim(integer_to_string(code)) &
00041          // '): ' // trim(warning_msg)
00042 
00043   end subroutine warn_msg
00044 
00045 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00046 
00047   !> Prints a warning message if condition_ok is false.
00048   subroutine warn_assert_msg(code, condition_ok, warning_msg)
00049 
00050     !> Status code to use.
00051     integer, intent(in) :: code
00052     !> Whether the assertion is ok.
00053     logical, intent(in) :: condition_ok
00054     !> Message to display.
00055     character(len=*), intent(in) :: warning_msg
00056 
00057     if (.not. condition_ok) then
00058        call warn_msg(code, warning_msg)
00059     end if
00060 
00061   end subroutine warn_assert_msg
00062 
00063 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00064 
00065   !> Errors unless condition_ok is true.
00066   subroutine assert_msg(code, condition_ok, error_msg)
00067 
00068     !> Status code to use if assertion fails.
00069     integer, intent(in) :: code
00070     !> Whether the assertion is ok.
00071     logical, intent(in) :: condition_ok
00072     !> Msg if assertion fails.
00073     character(len=*), intent(in) :: error_msg
00074 
00075     integer :: ierr
00076     character(len=100) :: code_str
00077 
00078     write(code_str,*) code
00079     code_str = adjustl(code_str)
00080     if (.not. condition_ok) then
00081        write(0,'(a)') 'ERROR (PartMC-' // trim(integer_to_string(code)) &
00082             // '): ' // trim(error_msg)
00083 #ifdef PMC_USE_MPI
00084        call mpi_abort(MPI_COMM_WORLD, code, ierr)
00085 #else
00086        stop 3
00087 #endif
00088     end if
00089 
00090   end subroutine assert_msg
00091 
00092 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00093 
00094   !> Errors unless condition_ok is true.
00095   subroutine assert(code, condition_ok)
00096 
00097     !> Status code to use if assertion fails.
00098     integer, intent(in) :: code
00099     !> Whether the assertion is ok.
00100     logical, intent(in) :: condition_ok
00101 
00102     if (.not. condition_ok) then
00103        ! much faster with gfortran 4.4.5 to do this extra test
00104        call assert_msg(code, condition_ok, 'assertion failed')
00105     end if
00106 
00107   end subroutine assert
00108 
00109 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00110 
00111   !> Error immediately.
00112   subroutine die(code)
00113 
00114     !> Status code to use if assertion fails.
00115     integer, intent(in) :: code
00116 
00117     call assert(code, .false.)
00118 
00119   end subroutine die
00120 
00121 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00122 
00123   !> Error immediately.
00124   subroutine die_msg(code, error_msg)
00125 
00126     !> Status code to use if assertion fails.
00127     integer, intent(in) :: code
00128     !> Msg if assertion fails.
00129     character(len=*), intent(in) :: error_msg
00130 
00131     call assert_msg(code, .false., error_msg)
00132 
00133   end subroutine die_msg
00134 
00135 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00136 
00137   !> Returns an available unit number. This should be freed by free_unit().
00138   integer function get_unit()
00139 
00140     integer i
00141     logical found_unit
00142 
00143     found_unit = .false.
00144     do i = 1,max_units
00145        if (.not. unit_used(i)) then
00146           found_unit = .true.
00147           exit
00148        end if
00149     end do
00150     if (.not. found_unit) then
00151        call die_msg(690355443, &
00152             'no more units available - need to free_unit()')
00153     end if
00154     unit_used(i) = .true.
00155     get_unit = i + unit_offset
00156 
00157   end function get_unit
00158 
00159 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00160 
00161   !> Frees a unit number returned by get_unit().
00162   subroutine free_unit(unit)
00163 
00164     integer, intent(in) :: unit
00165 
00166     unit_used(unit - unit_offset) = .false.
00167 
00168   end subroutine free_unit
00169 
00170 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00171 
00172   !> Convert volume (m^3) to radius (m).
00173   real(kind=dp) elemental function vol2rad(v)
00174 
00175     !> Volume (m^3).
00176     real(kind=dp), intent(in) :: v
00177     
00178     vol2rad = (v / (4d0 / 3d0 * const%pi))**(1d0/3d0)
00179     
00180   end function vol2rad
00181   
00182 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00183 
00184   !> Convert volume (m^3) to diameter (m).
00185   real(kind=dp) elemental function vol2diam(v)
00186 
00187     !> Volume (m^3).
00188     real(kind=dp), intent(in) :: v
00189     
00190     vol2diam = 2d0 * (v / (4d0 / 3d0 * const%pi))**(1d0/3d0)
00191     
00192   end function vol2diam
00193   
00194 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00195 
00196   !> Convert radius (m) to diameter (m).
00197   real(kind=dp) elemental function rad2diam(r)
00198 
00199     !> Radius (m).
00200     real(kind=dp), intent(in) :: r
00201     
00202     rad2diam = 2d0 * r
00203     
00204   end function rad2diam
00205   
00206 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00207 
00208   !> Convert radius (m) to volume (m^3).
00209   real(kind=dp) elemental function rad2vol(r)
00210 
00211     !> Radius (m).
00212     real(kind=dp), intent(in) :: r
00213     
00214     rad2vol = 4d0 / 3d0 * const%pi * r**3d0
00215     
00216   end function rad2vol
00217   
00218 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00219 
00220   !> Convert diameter (m) to radius (m).
00221   real(kind=dp) elemental function diam2rad(d)
00222 
00223     !> Diameter (m).
00224     real(kind=dp), intent(in) :: d
00225     
00226     diam2rad = d / 2d0
00227     
00228   end function diam2rad
00229   
00230 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00231 
00232   !> Convert diameter (m) to volume (m^3).
00233   real(kind=dp) elemental function diam2vol(d)
00234 
00235     !> Diameter (m).
00236     real(kind=dp), intent(in) :: d
00237     
00238     diam2vol = 4d0 / 3d0 * const%pi * (d / 2d0)**3d0
00239     
00240   end function diam2vol
00241   
00242 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00243 
00244   !> Tests whether two real numbers are almost equal using only a
00245   !> relative tolerance.
00246   logical function almost_equal(d1, d2)
00247     
00248     !> First number to compare.
00249     real(kind=dp), intent(in) :: d1
00250     !> Second number to compare.
00251     real(kind=dp), intent(in) :: d2
00252     
00253     !> Relative tolerance.
00254     real(kind=dp), parameter :: eps = 1d-8
00255     
00256     ! handle the 0.0 case
00257     if (d1 .eq. d2) then
00258        almost_equal = .true.
00259     else
00260        if (abs(d1 - d2) / (abs(d1) + abs(d2)) .lt. eps) then
00261           almost_equal = .true.
00262        else
00263           almost_equal = .false.
00264        end if
00265     end if
00266     
00267   end function almost_equal
00268   
00269 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00270 
00271   !> Tests whether two real numbers are almost equal using an
00272   !> absolute and relative tolerance.
00273   logical function almost_equal_abs(d1, d2, abs_tol)
00274     
00275     !> First number to compare.
00276     real(kind=dp), intent(in) :: d1
00277     !> Second number to compare.
00278     real(kind=dp), intent(in) :: d2
00279     !> Tolerance for when d1 equals d2.
00280     real(kind=dp), intent(in) :: abs_tol
00281     
00282     !> Relative tolerance.
00283     real(kind=dp), parameter :: eps = 1d-8
00284     
00285     ! handle the 0.0 case
00286     if (d1 .eq. d2) then
00287        almost_equal_abs = .true.
00288     else
00289        if ((abs(d1 - d2) .lt. abs_tol) .or. &
00290             (abs(d1 - d2) / (abs(d1) + abs(d2)) .lt. eps)) then
00291           almost_equal_abs = .true.
00292        else
00293           almost_equal_abs = .false.
00294        end if
00295     end if
00296     
00297   end function almost_equal_abs
00298   
00299 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00300 
00301   !> Computes whether an event is scheduled to take place.
00302   !!
00303   !! The events should occur ideally at times 0, interval, 2*interval,
00304   !! etc. The events are guaranteed to occur at least interval * (1 -
00305   !! tolerance) apart, and if at least interval time has passed then
00306   !! the next call is guaranteed to do the event. Otherwise the
00307   !! timestep is used to guess whether to do the event.
00308   subroutine check_event(time, timestep, interval, last_time, &
00309        do_event)
00310     
00311     !> Current time.
00312     real(kind=dp), intent(in) :: time
00313     !> Estimate of the time to the next call.
00314     real(kind=dp), intent(in) :: timestep
00315     !> How often the event should be done.
00316     real(kind=dp), intent(in) :: interval
00317     !> When the event was last done.
00318     real(kind=dp), intent(inout) :: last_time
00319     !> Whether the event should be done.
00320     logical, intent(out) :: do_event
00321     
00322     !> Fuzz for event occurance.
00323     real(kind=dp), parameter :: tolerance = 1d-6
00324     
00325     real(kind=dp) closest_interval_time
00326     
00327     ! if we are at time 0 then do the event unconditionally
00328     if (time .eq. 0d0) then
00329        do_event = .true.
00330     else
00331        ! if we are too close to the last time then don't do it
00332        if ((time - last_time) .lt. interval * (1d0 - tolerance)) then
00333           do_event = .false.
00334        else
00335           ! if it's been too long since the last time then do it
00336           if ((time - last_time) .ge. interval) then
00337              do_event = .true.
00338           else
00339              ! gray area -- if we are closer than we will be next
00340              ! time then do it
00341              closest_interval_time = anint(time / interval) * interval
00342              if (abs(time - closest_interval_time) &
00343                   .lt. abs(time + timestep - closest_interval_time)) &
00344                   then
00345                 do_event = .true.
00346              else
00347                 do_event = .false.
00348              end if
00349           end if
00350        end if
00351     end if
00352     
00353     if (do_event) then
00354        last_time = time
00355     end if
00356     
00357   end subroutine check_event
00358   
00359 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00360 
00361   !> Makes a linearly spaced array from min to max.
00362   subroutine linspace(min_x, max_x, x)
00363 
00364     !> Minimum array value.
00365     real(kind=dp), intent(in) :: min_x
00366     !> Maximum array value.
00367     real(kind=dp), intent(in) :: max_x
00368     !> Array.
00369     real(kind=dp), intent(out) :: x(:)
00370 
00371     integer :: i, n
00372     real(kind=dp) :: a
00373 
00374     n = size(x)
00375     do i = 2, (n - 1)
00376        a = real(i - 1, kind=dp) / real(n - 1, kind=dp)
00377        x(i) = (1d0 - a) * min_x + a * max_x
00378     end do
00379     if (n > 0) then
00380        ! make sure these values are exact
00381        x(1) = min_x
00382        x(n) = max_x
00383     end if
00384     
00385   end subroutine linspace
00386 
00387 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00388 
00389   !> Makes a logarithmically spaced array of length n from min to max.
00390   subroutine logspace(min_x, max_x, x)
00391 
00392     !> Minimum array value.
00393     real(kind=dp), intent(in) :: min_x
00394     !> Maximum array value.
00395     real(kind=dp), intent(in) :: max_x
00396     !> Array.
00397     real(kind=dp), intent(out) :: x(:)
00398 
00399     integer :: n
00400     real(kind=dp) :: log_x(size(x))
00401 
00402     n = size(x)
00403     call assert(548290438, min_x > 0d0)
00404     call assert(805259035, max_x > 0d0)
00405     call linspace(log(min_x), log(max_x), log_x)
00406     x = exp(log_x)
00407     if (n > 0) then
00408        ! make sure these values are exact
00409        x(1) = min_x
00410        x(n) = max_x
00411     end if
00412     
00413   end subroutine logspace
00414 
00415 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00416 
00417   !> Find the position of a real number in a 1D linear array.
00418   !!
00419   !! If xa is the array allocated by linspace(min_x, max_x, xa)
00420   !! then i = linspace_find(min_x, max_x, n, x) returns the index i
00421   !! satisfying xa(i) <= x < xa(i+1) for min_x <= x < max_x. If
00422   !! x >= max_x then i = n - 1.  If x < min_x then i = 1. Thus
00423   !! 1 <= i <= n - 1. Here n is the length of xa.
00424   !!
00425   !! This is equivalent to using find_1d() but much faster if the
00426   !! array is linear.
00427   integer function linspace_find(min_x, max_x, n, x)
00428 
00429     !> Minimum array value.
00430     real(kind=dp), intent(in) :: min_x
00431     !> Maximum array value.
00432     real(kind=dp), intent(in) :: max_x
00433     !> Number of entries.
00434     integer, intent(in) :: n
00435     !> Value.
00436     real(kind=dp), intent(in) :: x
00437 
00438     linspace_find = floor((x - min_x) / (max_x - min_x) &
00439          * real(n - 1, kind=dp)) + 1
00440     linspace_find = min(linspace_find, n - 1)
00441     linspace_find = max(linspace_find, 1)
00442     
00443   end function linspace_find
00444 
00445 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00446 
00447   !> Find the position of a real number in a 1D logarithmic array.
00448   !!
00449   !! If xa is the array allocated by logspace(min_x, max_x, n, xa) then
00450   !! i = logspace_find(min_x, max_x, n, x) returns the index i satisfying
00451   !! xa(i) <= x < xa(i+1) for min_x <= x < max_x. If x >= max_x then i = n.
00452   !! If x < min_x then i = 1. Thus 1 <= i <= n.
00453   !!
00454   !! This is equivalent to using find_1d() but much faster if the
00455   !! array is known to be logarithmic.
00456   integer function logspace_find(min_x, max_x, n, x)
00457 
00458     !> Minimum array value.
00459     real(kind=dp), intent(in) :: min_x
00460     !> Maximum array value.
00461     real(kind=dp), intent(in) :: max_x
00462     !> Number of entries.
00463     integer, intent(in) :: n
00464     !> Value.
00465     real(kind=dp), intent(in) :: x
00466 
00467     logspace_find = linspace_find(log(min_x), log(max_x), n, log(x))
00468     
00469   end function logspace_find
00470 
00471 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00472 
00473   !> Find the position of a real number in an arbitrary 1D array.
00474   !!
00475   !! Takes an array of x_vals, and a single x value, and returns the
00476   !! position p such that x_vals(p) <= x < x_vals(p+1). If p == 0 then
00477   !! x < x_vals(1) and if p == n then x_vals(n) <= x. x_vals must be
00478   !! sorted in increasing order.
00479   !!
00480   !! If the array is known to be linearly or logarithmically spaced
00481   !! then linspace_find() or logspace_find() will be much faster.
00482   integer function find_1d(n, x_vals, x)
00483 
00484     !> Number of values.
00485     integer, intent(in) :: n
00486     !> X value array, must be sorted.
00487     real(kind=dp), intent(in) :: x_vals(n)
00488     !> Value to interpolate at.
00489     real(kind=dp), intent(in) :: x
00490 
00491     integer p
00492 
00493     if (n == 0) then
00494        find_1d = 0
00495        return
00496     end if
00497     p = 1
00498     do while (x >= x_vals(p))
00499        p = p + 1
00500        if (p > n) then
00501           exit
00502        end if
00503     end do
00504     p = p - 1
00505     find_1d = p
00506 
00507   end function find_1d
00508 
00509 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00510 
00511   !> 1D linear interpolation.
00512   !!
00513   !! Takes an array of x and y, and a single x value, and returns the
00514   !! corresponding y using linear interpolation. x_vals must be
00515   !! sorted.
00516   real(kind=dp) function interp_1d(x_vals, y_vals, x)
00517 
00518     !> X value array, must be sorted.
00519     real(kind=dp), intent(in) :: x_vals(:)
00520     !> Y value array.
00521     real(kind=dp), intent(in) :: y_vals(size(x_vals))
00522     !> Value to interpolate at.
00523     real(kind=dp), intent(in) :: x
00524 
00525     integer :: n, p
00526     real(kind=dp) :: y, alpha
00527 
00528     n = size(x_vals)
00529     p = find_1d(n, x_vals, x)
00530     if (p == 0) then
00531        y = y_vals(1)
00532     elseif (p == n) then
00533        y = y_vals(n)
00534     else
00535        alpha = (x - x_vals(p)) / (x_vals(p+1) - x_vals(p))
00536        y = (1d0 - alpha) * y_vals(p) + alpha * y_vals(p+1)
00537     end if
00538     interp_1d = y
00539 
00540   end function interp_1d
00541 
00542 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00543 
00544   !> Linear interpolation over discrete indices.
00545   !!
00546   !! Takes real values \c x_1 and \c x_n and integers \c n and \c i
00547   !! and returns the linear interpolation so that \c x_1 is returned
00548   !! when \c i = 1 and \c x_n is returned when \c i = \c n.
00549   real(kind=dp) function interp_linear_disc(x_1, x_n, n, i)
00550 
00551     !> Value of \c x when \c i = 1.
00552     real(kind=dp), intent(in) :: x_1
00553     !> Value of \c x when \c i = n.
00554     real(kind=dp), intent(in) :: x_n
00555     !> Number of points to interpolate over.
00556     integer, intent(in) :: n
00557     !> Index to interpolate at.
00558     integer, intent(in) :: i
00559 
00560     real(kind=dp) :: alpha
00561 
00562     alpha = real(i - 1, kind=dp) / real(n - 1, kind=dp)
00563     interp_linear_disc = (1d0 - alpha) * x_1 + alpha * x_n
00564 
00565   end function interp_linear_disc
00566 
00567 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00568 
00569   !> Convert a string to an integer.
00570   integer function string_to_integer(string)
00571 
00572     !> String to convert.
00573     character(len=*), intent(in) :: string
00574     
00575     integer :: val
00576     integer :: ios
00577     character(len=len(string)+300) :: error_msg
00578 
00579     read(string, '(i20)', iostat=ios) val
00580     call assert_msg(895881873, ios == 0, &
00581          'error converting "' // trim(string) &
00582          // '" to integer: IOSTAT = ' // trim(integer_to_string(ios)))
00583     string_to_integer = val
00584 
00585   end function string_to_integer
00586 
00587 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00588 
00589   !> Convert a string to a real.
00590   real(kind=dp) function string_to_real(string)
00591 
00592     !> String to convert.
00593     character(len=*), intent(in) :: string
00594     
00595     real(kind=dp) :: val
00596     integer :: ios
00597     character(len=len(string)+300) :: error_msg
00598 
00599     read(string, '(f20.0)', iostat=ios) val
00600     call assert_msg(727430097, ios == 0, &
00601          'error converting "' // trim(string) &
00602          // '" to real: IOSTAT = ' // trim(integer_to_string(ios)))
00603     string_to_real = val
00604 
00605   end function string_to_real
00606 
00607 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00608 
00609   !> Convert a string to a logical.
00610   logical function string_to_logical(string)
00611 
00612     !> String to convert.
00613     character(len=*), intent(in) :: string
00614     
00615     logical :: val
00616     integer :: ios
00617     character(len=len(string)+300) :: error_msg
00618 
00619     val = .false.
00620     if ((trim(string) == 'yes') &
00621          .or. (trim(string) == 'y') &
00622          .or. (trim(string) == 'true') &
00623          .or. (trim(string) == 't') &
00624          .or. (trim(string) == '1')) then
00625        val = .true.
00626     elseif ((trim(string) == 'no') &
00627          .or. (trim(string) == 'n') &
00628          .or. (trim(string) == 'false') &
00629          .or. (trim(string) == 'f') &
00630          .or. (trim(string) == '0')) then
00631        val = .false.
00632     else
00633        call assert_msg(985010153, ios == 0, &
00634             'error converting "' // trim(string) &
00635             // '" to logical: IOSTAT = ' // trim(integer_to_string(ios)))
00636     end if
00637     string_to_logical = val
00638 
00639   end function string_to_logical
00640 
00641 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00642 
00643   !> Convert an integer to a string format.
00644   character(len=PMC_UTIL_CONVERT_STRING_LEN) function integer_to_string(val)
00645 
00646     !> Value to convert.
00647     integer, intent(in) :: val
00648 
00649     character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
00650     
00651     ret_val = ""
00652     write(ret_val, '(i30)') val
00653     integer_to_string = adjustl(ret_val)
00654 
00655   end function integer_to_string
00656 
00657 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00658 
00659   !> Convert a real to a string format.
00660   character(len=PMC_UTIL_CONVERT_STRING_LEN) function real_to_string(val)
00661 
00662     !> Value to convert.
00663     real(kind=dp), intent(in) :: val
00664 
00665     character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
00666     
00667     ret_val = ""
00668     write(ret_val, '(g30.20)') val
00669     real_to_string = adjustl(ret_val)
00670 
00671   end function real_to_string
00672 
00673 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00674 
00675   !> Convert a logical to a string format.
00676   character(len=PMC_UTIL_CONVERT_STRING_LEN) function logical_to_string(val)
00677 
00678     !> Value to convert.
00679     logical, intent(in) :: val
00680 
00681     character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
00682     
00683     ret_val = ""
00684     if (val) then
00685        ret_val = "TRUE"
00686     else
00687        ret_val = "FALSE"
00688     end if
00689     logical_to_string = ret_val
00690 
00691   end function logical_to_string
00692 
00693 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00694 
00695   !> Convert a complex to a string format.
00696   character(len=PMC_UTIL_CONVERT_STRING_LEN) function complex_to_string(val)
00697 
00698     !> Value to convert.
00699     complex(kind=dc), intent(in) :: val
00700 
00701     character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
00702 
00703     ret_val = ""
00704     ret_val = "(" // trim(real_to_string(real(val))) 
00705          // ", " // trim(real_to_string(aimag(val))) // ")"
00706     complex_to_string = adjustl(ret_val)
00707 
00708   end function complex_to_string
00709 
00710 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00711 
00712   !> Convert an integer to a string format of maximum length.
00713   character(len=PMC_UTIL_CONVERT_STRING_LEN) 
00714        function integer_to_string_max_len(val, max_len)
00715 
00716     !> Value to convert.
00717     integer, intent(in) :: val
00718     !> Maximum length of resulting string.
00719     integer, intent(in) :: max_len
00720 
00721     character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
00722 
00723     ret_val = integer_to_string(val)
00724     if (len_trim(ret_val) > max_len) then
00725        ret_val = real_to_string_max_len(real(val, kind=dp), max_len)
00726     end if
00727     integer_to_string_max_len = adjustl(ret_val)
00728 
00729   end function integer_to_string_max_len
00730 
00731 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00732 
00733   !> Convert a real to a string format of maximum length.
00734   character(len=PMC_UTIL_CONVERT_STRING_LEN) 
00735        function real_to_string_max_len(val, max_len)
00736 
00737     !> Value to convert.
00738     real(kind=dp), intent(in) :: val
00739     !> Maximum length of resulting string.
00740     integer, intent(in) :: max_len
00741 
00742     character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val, exp_str, frac_str
00743     integer :: exp_val, exp_len, frac_len, use_frac_len, min_frac_len, i
00744     real(kind=dp) :: frac_val
00745 
00746     ret_val = ""
00747     if (val == 0d0) then
00748        if (max_len >= 3) then
00749           ret_val = "0e0"
00750        else
00751           do i = 1,max_len
00752              ret_val(i:i) = "*"
00753           end do
00754        end if
00755        real_to_string_max_len = adjustl(ret_val)
00756        return
00757     end if
00758 
00759     exp_val = floor(log10(abs(val)))
00760     frac_val = val / 10d0**exp_val
00761     exp_str = integer_to_string(exp_val)
00762     frac_str = real_to_string(frac_val)
00763 
00764     exp_len = len_trim(exp_str)
00765     frac_len = len_trim(frac_str)
00766     use_frac_len = max_len - 1 - exp_len
00767     if (use_frac_len > frac_len) then
00768        use_frac_len = frac_len
00769     end if
00770     if (val < 0d0) then
00771        min_frac_len = 2
00772     else
00773        min_frac_len = 1
00774     end if
00775     if (use_frac_len < min_frac_len) then
00776        do i = 1,max_len
00777           ret_val(i:i) = "*"
00778        end do
00779     else
00780        ret_val = frac_str(1:use_frac_len) // "e" // trim(exp_str)
00781     end if
00782     real_to_string_max_len = adjustl(ret_val)
00783 
00784   end function real_to_string_max_len
00785 
00786 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00787 
00788   !> Convert a time to a string format of maximum length.
00789   character(len=PMC_UTIL_CONVERT_STRING_LEN) 
00790        function time_to_string_max_len(time, max_len)
00791 
00792     !> Time to convert (s).
00793     real(kind=dp), intent(in) :: time
00794     !> Maximum length of resulting string.
00795     integer, intent(in) :: max_len
00796 
00797     integer, dimension(4), parameter :: scale  = (/   1,  60,  60,  24 /)
00798     character, dimension(4), parameter :: unit = (/ "s", "m", "h", "d" /)
00799 
00800     character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
00801     integer :: i
00802     logical :: len_ok
00803     real(kind=dp) :: scaled_time
00804 
00805     scaled_time = time
00806     len_ok = .false.
00807     do i = 1,4
00808        scaled_time = scaled_time / real(scale(i), kind=dp)
00809        ret_val = trim(integer_to_string(nint(scaled_time))) // unit(i)
00810        if (len_trim(ret_val) <= max_len) then
00811           len_ok = .true.
00812           exit
00813        end if
00814     end do
00815     if (.not. len_ok) then
00816        ret_val = ""
00817        do i = 1,max_len
00818           ret_val(i:i) = "*"
00819        end do
00820     end if
00821     time_to_string_max_len = adjustl(ret_val)
00822 
00823   end function time_to_string_max_len
00824 
00825 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00826 
00827   !> Convert a real-valued vector into an integer-valued vector.
00828   !!
00829   !! Use n_samp samples. Returns discrete vector whose relative entry
00830   !! sizes are \f$ \ell_1 \f$ closest to the continuous vector.
00831   subroutine vec_cts_to_disc(n, vec_cts, n_samp, vec_disc)
00832     
00833     !> Number of entries in vectors.
00834     integer, intent(in) :: n
00835     !> Continuous vector.
00836     real(kind=dp), intent(in) :: vec_cts(n)
00837     !> Number of discrete samples to use.
00838     integer, intent(in) :: n_samp
00839     !> Discretized vector.
00840     integer, intent(out) :: vec_disc(n)
00841     
00842     integer :: k(1)
00843     real(kind=dp) :: vec_tot
00844     
00845     vec_tot = sum(vec_cts)
00846     
00847     ! assign a best guess for each bin independently
00848     vec_disc = nint(vec_cts / vec_tot * real(n_samp, kind=dp))
00849     
00850     ! if we have too few particles then add more
00851     do while (sum(vec_disc) < n_samp)
00852        k = minloc(abs(real(vec_disc + 1, kind=dp) - vec_cts) 
00853             - abs(real(vec_disc, kind=dp) - vec_cts))
00854        vec_disc(k) = vec_disc(k) + 1
00855     end do
00856     
00857     ! if we have too many particles then remove some
00858     do while (sum(vec_disc) > n_samp)
00859        k = minloc(abs(real(vec_disc - 1, kind=dp) - vec_cts) 
00860             - abs(real(vec_disc, kind=dp) - vec_cts))
00861        vec_disc(k) = vec_disc(k) - 1
00862     end do
00863     
00864     call assert_msg(323412496, sum(vec_disc) == n_samp, &
00865          'generated incorrect number of samples')
00866     
00867   end subroutine vec_cts_to_disc
00868   
00869 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00870 
00871   !> Computes the average of an array of integer numbers.
00872   subroutine average_integer(int_vec, int_avg)
00873     
00874     !> Array of integer numbers.
00875     integer, intent(in) :: int_vec(:)
00876     !> Average of int_vec.
00877     integer, intent(out) :: int_avg
00878     
00879     int_avg = sum(int_vec) / size(int_vec)
00880     
00881   end subroutine average_integer
00882   
00883 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00884 
00885   !> Computes the average of an array of real numbers.
00886   subroutine average_real(real_vec, real_avg)
00887     
00888     !> Array of real numbers.
00889     real(kind=dp), intent(in) :: real_vec(:)
00890     !> Average of real_vec.
00891     real(kind=dp), intent(out) :: real_avg
00892     
00893     real_avg = sum(real_vec) / real(size(real_vec), kind=dp)
00894     
00895   end subroutine average_real
00896   
00897 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00898 
00899   !> Strip the extension to find the basename.
00900   !!
00901   !! The filename is assumed to be of the form basename.extension,
00902   !! where extension contains no periods. If there is no period in
00903   !! filename then basename is the whole filename.
00904   subroutine get_basename(filename, basename)
00905 
00906     !> Full filename.
00907     character(len=*), intent(in) :: filename
00908     !> Basename.
00909     character(len=*), intent(out) :: basename
00910 
00911     integer :: i
00912     logical :: found_period
00913 
00914     basename = filename
00915     i = len_trim(basename)
00916     found_period = .false.
00917     do while ((i > 0) .and. (.not. found_period))
00918        ! Fortran .and. does not short-circuit, so we can't do the
00919        ! obvious do while ((i > 0) .and. (basename(i:i) /= ".")),
00920        ! instead we have to use this hack with the found_period
00921        ! logical variable.
00922        if (basename(i:i) == ".") then
00923           found_period = .true.
00924        else
00925           i = i - 1
00926        end if
00927     end do
00928     if (i > 0) then
00929        basename(i:) = ""
00930     end if
00931 
00932   end subroutine get_basename
00933 
00934 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00935 
00936   !> Current date and time in ISO 8601 format.
00937   subroutine iso8601_date_and_time(date_time)
00938 
00939     !> Result string.
00940     character(len=*), intent(out) :: date_time
00941 
00942     character(len=10) :: date, time, zone
00943 
00944     call assert_msg(893219839, len(date_time) >= 29, &
00945          "date_time string must have length at least 29")
00946     call date_and_time(date, time, zone)
00947     date_time = ""
00948     write(date_time, '(14a)') date(1:4), "-", date(5:6), "-", &
00949          date(7:8), "T", time(1:2), ":", time(3:4), ":", &
00950          time(5:10), zone(1:3), ":", zone(4:5)
00951     
00952   end subroutine iso8601_date_and_time
00953 
00954 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00955 
00956   !> Convert degrees to radians.
00957   real(kind=dp) function deg2rad(deg)
00958 
00959     !> Input degrees.
00960     real(kind=dp), intent(in) :: deg
00961 
00962     deg2rad = deg / 180d0 * const%pi
00963 
00964   end function deg2rad
00965 
00966 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00967 
00968   !> Convert radians to degrees.
00969   real(kind=dp) function rad2deg(rad)
00970 
00971     !> Input radians.
00972     real(kind=dp), intent(in) :: rad
00973 
00974     rad2deg = rad / const%pi * 180d0
00975 
00976   end function rad2deg
00977   
00978 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00979   
00980 #ifdef DEFINE_LOCAL_COMMAND_ARGUMENT
00981   !> Hack for compilers that don't support the Fortran2003 standard
00982   !> command_argument_count() function.
00983   integer function command_argument_count()
00984 
00985     command_argument_count = iargc()
00986 
00987   end function command_argument_count
00988 #endif
00989 
00990 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00991 
00992 #ifdef DEFINE_LOCAL_COMMAND_ARGUMENT
00993   !> Hack for compilers that don't support the Fortran2003 standard
00994   !> get_command_argument() subroutine.
00995   subroutine get_command_argument(i, arg)
00996     integer, intent(in) :: i
00997     character(len=*), intent(out) :: arg
00998 
00999     call getarg(i, arg)
01000 
01001   end subroutine get_command_argument
01002 #endif
01003   
01004 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
01005 
01006 end module pmc_util