PartMC
2.2.0
|
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 #ifdef PMC_USE_C_SORT 00016 use iso_c_binding 00017 #endif 00018 00019 !> Maximum number of IO units usable simultaneously. 00020 integer, parameter :: max_units = 200 00021 !> Minimum unit number to allocate. 00022 integer, parameter :: unit_offset = 10 00023 !> Table of unit numbers storing allocation status. 00024 logical, save :: unit_used(max_units) = .false. 00025 00026 !> Length of string for converting numbers. 00027 integer, parameter :: PMC_UTIL_CONVERT_STRING_LEN = 100 00028 !> Maximum length of filenames. 00029 integer, parameter :: PMC_MAX_FILENAME_LEN = 300 00030 00031 contains 00032 00033 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00034 00035 !> Prints a warning message. 00036 subroutine warn_msg(code, warning_msg) 00037 00038 !> Status code to use. 00039 integer, intent(in) :: code 00040 !> Message to display. 00041 character(len=*), intent(in) :: warning_msg 00042 00043 write(0,'(a)') 'WARNING (PartMC-' // trim(integer_to_string(code)) & 00044 // '): ' // trim(warning_msg) 00045 00046 end subroutine warn_msg 00047 00048 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00049 00050 !> Prints a warning message if condition_ok is false. 00051 subroutine warn_assert_msg(code, condition_ok, warning_msg) 00052 00053 !> Status code to use. 00054 integer, intent(in) :: code 00055 !> Whether the assertion is ok. 00056 logical, intent(in) :: condition_ok 00057 !> Message to display. 00058 character(len=*), intent(in) :: warning_msg 00059 00060 if (.not. condition_ok) then 00061 call warn_msg(code, warning_msg) 00062 end if 00063 00064 end subroutine warn_assert_msg 00065 00066 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00067 00068 !> Errors unless condition_ok is true. 00069 subroutine assert_msg(code, condition_ok, error_msg) 00070 00071 !> Status code to use if assertion fails. 00072 integer, intent(in) :: code 00073 !> Whether the assertion is ok. 00074 logical, intent(in) :: condition_ok 00075 !> Msg if assertion fails. 00076 character(len=*), intent(in) :: error_msg 00077 00078 integer :: ierr 00079 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 ! FIXME: is it still faster now that assert_msg doesn't 00105 ! unconditionally construct a code_str? 00106 call assert_msg(code, condition_ok, 'assertion failed') 00107 end if 00108 00109 end subroutine assert 00110 00111 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00112 00113 !> Error immediately. 00114 subroutine die(code) 00115 00116 !> Status code to use if assertion fails. 00117 integer, intent(in) :: code 00118 00119 call assert(code, .false.) 00120 00121 end subroutine die 00122 00123 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00124 00125 !> Error immediately. 00126 subroutine die_msg(code, error_msg) 00127 00128 !> Status code to use if assertion fails. 00129 integer, intent(in) :: code 00130 !> Msg if assertion fails. 00131 character(len=*), intent(in) :: error_msg 00132 00133 call assert_msg(code, .false., error_msg) 00134 00135 end subroutine die_msg 00136 00137 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00138 00139 !> Returns an available unit number. This should be freed by free_unit(). 00140 integer function get_unit() 00141 00142 integer i 00143 logical found_unit 00144 00145 found_unit = .false. 00146 do i = 1,max_units 00147 if (.not. unit_used(i)) then 00148 found_unit = .true. 00149 exit 00150 end if 00151 end do 00152 if (.not. found_unit) then 00153 call die_msg(690355443, & 00154 'no more units available - need to free_unit()') 00155 end if 00156 unit_used(i) = .true. 00157 get_unit = i + unit_offset 00158 00159 end function get_unit 00160 00161 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00162 00163 !> Frees a unit number returned by get_unit(). 00164 subroutine free_unit(unit) 00165 00166 integer, intent(in) :: unit 00167 00168 unit_used(unit - unit_offset) = .false. 00169 00170 end subroutine free_unit 00171 00172 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00173 00174 !> Open a file for reading with an automatically assigned unit and 00175 !> test that it succeeds. The file should be closed with 00176 !> close_file(). 00177 subroutine open_file_read(filename, unit) 00178 00179 !> Filename to open. 00180 character(len=*), intent(in) :: filename 00181 !> Unit assigned to file. 00182 integer, intent(out) :: unit 00183 00184 integer :: ios 00185 00186 unit = get_unit() 00187 open(unit=unit, file=filename, status='old', action='read', iostat=ios) 00188 call assert_msg(609624199, ios == 0, 'unable to open file ' & 00189 // trim(filename) // ' for reading: ' // integer_to_string(ios)) 00190 00191 end subroutine open_file_read 00192 00193 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00194 00195 !> Open a file for writing with an automatically assigned unit and 00196 !> test that it succeeds. The file should be closed with 00197 !> close_file(). 00198 subroutine open_file_write(filename, unit) 00199 00200 !> Filename to open. 00201 character(len=*), intent(in) :: filename 00202 !> Unit assigned to file. 00203 integer, intent(out) :: unit 00204 00205 integer :: ios 00206 00207 unit = get_unit() 00208 open(unit=unit, file=filename, status='replace', action='write', & 00209 iostat=ios) 00210 call assert_msg(609624199, ios == 0, 'unable to open file ' & 00211 // trim(filename) // ' for writing: ' // integer_to_string(ios)) 00212 00213 end subroutine open_file_write 00214 00215 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00216 00217 !> Close a file and de-assign the unit. 00218 subroutine close_file(unit) 00219 00220 !> Unit to close. 00221 integer, intent(in) :: unit 00222 00223 close(unit) 00224 call free_unit(unit) 00225 00226 end subroutine close_file 00227 00228 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00229 00230 !> Convert volume (m^3) to radius (m). 00231 real(kind=dp) elemental function vol2rad(v) 00232 00233 !> Volume (m^3). 00234 real(kind=dp), intent(in) :: v 00235 00236 vol2rad = (v / (4d0 / 3d0 * const%pi))**(1d0/3d0) 00237 00238 end function vol2rad 00239 00240 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00241 00242 !> Convert volume (m^3) to diameter (m). 00243 real(kind=dp) elemental function vol2diam(v) 00244 00245 !> Volume (m^3). 00246 real(kind=dp), intent(in) :: v 00247 00248 vol2diam = 2d0 * (v / (4d0 / 3d0 * const%pi))**(1d0/3d0) 00249 00250 end function vol2diam 00251 00252 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00253 00254 !> Convert radius (m) to diameter (m). 00255 real(kind=dp) elemental function rad2diam(r) 00256 00257 !> Radius (m). 00258 real(kind=dp), intent(in) :: r 00259 00260 rad2diam = 2d0 * r 00261 00262 end function rad2diam 00263 00264 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00265 00266 !> Convert radius (m) to volume (m^3). 00267 real(kind=dp) elemental function rad2vol(r) 00268 00269 !> Radius (m). 00270 real(kind=dp), intent(in) :: r 00271 00272 rad2vol = 4d0 / 3d0 * const%pi * r**3d0 00273 00274 end function rad2vol 00275 00276 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00277 00278 !> Convert diameter (m) to radius (m). 00279 real(kind=dp) elemental function diam2rad(d) 00280 00281 !> Diameter (m). 00282 real(kind=dp), intent(in) :: d 00283 00284 diam2rad = d / 2d0 00285 00286 end function diam2rad 00287 00288 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00289 00290 !> Convert diameter (m) to volume (m^3). 00291 real(kind=dp) elemental function diam2vol(d) 00292 00293 !> Diameter (m). 00294 real(kind=dp), intent(in) :: d 00295 00296 diam2vol = 4d0 / 3d0 * const%pi * (d / 2d0)**3d0 00297 00298 end function diam2vol 00299 00300 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00301 00302 !> Tests whether two real numbers are almost equal using only a 00303 !> relative tolerance. 00304 logical function almost_equal(d1, d2) 00305 00306 !> First number to compare. 00307 real(kind=dp), intent(in) :: d1 00308 !> Second number to compare. 00309 real(kind=dp), intent(in) :: d2 00310 00311 !> Relative tolerance. 00312 real(kind=dp), parameter :: eps = 1d-8 00313 00314 ! handle the 0.0 case 00315 if (d1 .eq. d2) then 00316 almost_equal = .true. 00317 else 00318 if (abs(d1 - d2) / (abs(d1) + abs(d2)) .lt. eps) then 00319 almost_equal = .true. 00320 else 00321 almost_equal = .false. 00322 end if 00323 end if 00324 00325 end function almost_equal 00326 00327 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00328 00329 !> Tests whether two real numbers are almost equal using an 00330 !> absolute and relative tolerance. 00331 logical function almost_equal_abs(d1, d2, abs_tol) 00332 00333 !> First number to compare. 00334 real(kind=dp), intent(in) :: d1 00335 !> Second number to compare. 00336 real(kind=dp), intent(in) :: d2 00337 !> Tolerance for when d1 equals d2. 00338 real(kind=dp), intent(in) :: abs_tol 00339 00340 !> Relative tolerance. 00341 real(kind=dp), parameter :: eps = 1d-8 00342 00343 ! handle the 0.0 case 00344 if (d1 .eq. d2) then 00345 almost_equal_abs = .true. 00346 else 00347 if ((abs(d1 - d2) .lt. abs_tol) .or. & 00348 (abs(d1 - d2) / (abs(d1) + abs(d2)) .lt. eps)) then 00349 almost_equal_abs = .true. 00350 else 00351 almost_equal_abs = .false. 00352 end if 00353 end if 00354 00355 end function almost_equal_abs 00356 00357 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00358 00359 !> Check that the first time interval is close to an integer 00360 !> multiple of the second, and warn if it is not. 00361 subroutine check_time_multiple(first_name, first_time, & 00362 second_name, second_time) 00363 00364 !> Name of the first time variable. 00365 character(len=*), intent(in) :: first_name 00366 !> First time variable (s). 00367 real(kind=dp), intent(in) :: first_time 00368 !> Name of the second time variable. 00369 character(len=*), intent(in) :: second_name 00370 !> Second time variable (s). 00371 real(kind=dp), intent(in) :: second_time 00372 00373 if (abs(mod(first_time, second_time) / second_time) > 1d-6) then 00374 call warn_msg(952299377, trim(first_name) & 00375 // " is not an integer multiple of " // trim(second_name)) 00376 end if 00377 00378 end subroutine check_time_multiple 00379 00380 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00381 00382 !> Computes whether an event is scheduled to take place. 00383 !! 00384 !! The events should occur ideally at times 0, interval, 2*interval, 00385 !! etc. The events are guaranteed to occur at least interval * (1 - 00386 !! tolerance) apart, and if at least interval time has passed then 00387 !! the next call is guaranteed to do the event. Otherwise the 00388 !! timestep is used to guess whether to do the event. 00389 subroutine check_event(time, timestep, interval, last_time, & 00390 do_event) 00391 00392 !> Current time. 00393 real(kind=dp), intent(in) :: time 00394 !> Estimate of the time to the next call. 00395 real(kind=dp), intent(in) :: timestep 00396 !> How often the event should be done. 00397 real(kind=dp), intent(in) :: interval 00398 !> When the event was last done. 00399 real(kind=dp), intent(inout) :: last_time 00400 !> Whether the event should be done. 00401 logical, intent(out) :: do_event 00402 00403 !> Fuzz for event occurance. 00404 real(kind=dp), parameter :: tolerance = 1d-6 00405 00406 real(kind=dp) closest_interval_time 00407 00408 ! if we are at time 0 then do the event unconditionally 00409 if (time .eq. 0d0) then 00410 do_event = .true. 00411 else 00412 ! if we are too close to the last time then don't do it 00413 if ((time - last_time) .lt. interval * (1d0 - tolerance)) then 00414 do_event = .false. 00415 else 00416 ! if it's been too long since the last time then do it 00417 if ((time - last_time) .ge. interval) then 00418 do_event = .true. 00419 else 00420 ! gray area -- if we are closer than we will be next 00421 ! time then do it 00422 closest_interval_time = anint(time / interval) * interval 00423 if (abs(time - closest_interval_time) & 00424 .lt. abs(time + timestep - closest_interval_time)) & 00425 then 00426 do_event = .true. 00427 else 00428 do_event = .false. 00429 end if 00430 end if 00431 end if 00432 end if 00433 00434 if (do_event) then 00435 last_time = time 00436 end if 00437 00438 end subroutine check_event 00439 00440 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00441 00442 !> Makes a linearly spaced array from min to max. 00443 subroutine linspace(min_x, max_x, x) 00444 00445 !> Minimum array value. 00446 real(kind=dp), intent(in) :: min_x 00447 !> Maximum array value. 00448 real(kind=dp), intent(in) :: max_x 00449 !> Array. 00450 real(kind=dp), intent(out) :: x(:) 00451 00452 integer :: i, n 00453 real(kind=dp) :: a 00454 00455 n = size(x) 00456 do i = 2, (n - 1) 00457 a = real(i - 1, kind=dp) / real(n - 1, kind=dp) 00458 x(i) = (1d0 - a) * min_x + a * max_x 00459 end do 00460 if (n > 0) then 00461 ! make sure these values are exact 00462 x(1) = min_x 00463 x(n) = max_x 00464 end if 00465 00466 end subroutine linspace 00467 00468 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00469 00470 !> Makes a logarithmically spaced array of length n from min to max. 00471 subroutine logspace(min_x, max_x, x) 00472 00473 !> Minimum array value. 00474 real(kind=dp), intent(in) :: min_x 00475 !> Maximum array value. 00476 real(kind=dp), intent(in) :: max_x 00477 !> Array. 00478 real(kind=dp), intent(out) :: x(:) 00479 00480 integer :: n 00481 real(kind=dp) :: log_x(size(x)) 00482 00483 n = size(x) 00484 call assert(585043176, n > 0) 00485 call assert(548290438, min_x > 0d0) 00486 call assert(805259035, max_x > 0d0) 00487 call linspace(log(min_x), log(max_x), log_x) 00488 x = exp(log_x) 00489 if (n > 0) then 00490 ! make sure these values are exact 00491 x(1) = min_x 00492 x(n) = max_x 00493 end if 00494 00495 end subroutine logspace 00496 00497 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00498 00499 !> Find the position of a real number in a 1D linear array. 00500 !! 00501 !! If xa is the array allocated by linspace(min_x, max_x, xa) then i 00502 !! = linspace_find(min_x, max_x, n, x) returns the index i 00503 !! satisfying xa(i) <= x < xa(i+1) for min_x <= x < max_x. If x >= 00504 !! max_x then i = n. If x < min_x then i = 0. Thus 0 <= i <= 00505 !! n. Here n is the length of xa. 00506 !! 00507 !! This is equivalent to using find_1d() but much faster if the 00508 !! array is linear. 00509 integer function linspace_find(min_x, max_x, n, x) 00510 00511 !> Minimum array value. 00512 real(kind=dp), intent(in) :: min_x 00513 !> Maximum array value. 00514 real(kind=dp), intent(in) :: max_x 00515 !> Number of entries. 00516 integer, intent(in) :: n 00517 !> Value. 00518 real(kind=dp), intent(in) :: x 00519 00520 linspace_find = floor((x - min_x) / (max_x - min_x) & 00521 * real(n - 1, kind=dp)) + 1 00522 linspace_find = min(linspace_find, n) 00523 linspace_find = max(linspace_find, 0) 00524 00525 end function linspace_find 00526 00527 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00528 00529 !> Find the position of a real number in a 1D logarithmic array. 00530 !! 00531 !! If xa is the array allocated by logspace(min_x, max_x, xa) then i 00532 !! = logspace_find(min_x, max_x, n, x) returns the index i 00533 !! satisfying xa(i) <= x < xa(i+1) for min_x <= x < max_x. If x >= 00534 !! max_x then i = n. If x < min_x then i = 0. Thus 0 <= i <= 00535 !! n. Here n is the length of xa. 00536 !! 00537 !! This is equivalent to using find_1d() but much faster if the 00538 !! array is logarithmic. 00539 integer function logspace_find(min_x, max_x, n, x) 00540 00541 !> Minimum array value. 00542 real(kind=dp), intent(in) :: min_x 00543 !> Maximum array value. 00544 real(kind=dp), intent(in) :: max_x 00545 !> Number of entries. 00546 integer, intent(in) :: n 00547 !> Value. 00548 real(kind=dp), intent(in) :: x 00549 00550 logspace_find = linspace_find(log(min_x), log(max_x), n, log(x)) 00551 00552 end function logspace_find 00553 00554 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00555 00556 !> Find the position of a real number in an arbitrary 1D array. 00557 !! 00558 !! Takes an array of x_vals, and a single x value, and returns the 00559 !! position p such that x_vals(p) <= x < x_vals(p+1). If p == 0 then 00560 !! x < x_vals(1) and if p == n then x_vals(n) <= x. x_vals must be 00561 !! sorted in increasing order. 00562 !! 00563 !! If the array is known to be linearly or logarithmically spaced 00564 !! then linspace_find() or logspace_find() will be much faster. 00565 integer function find_1d(n, x_vals, x) 00566 00567 !> Number of values. 00568 integer, intent(in) :: n 00569 !> X value array, must be sorted. 00570 real(kind=dp), intent(in) :: x_vals(n) 00571 !> Value to interpolate at. 00572 real(kind=dp), intent(in) :: x 00573 00574 integer p 00575 00576 if (n == 0) then 00577 find_1d = 0 00578 return 00579 end if 00580 p = 1 00581 do while (x >= x_vals(p)) 00582 p = p + 1 00583 if (p > n) then 00584 exit 00585 end if 00586 end do 00587 p = p - 1 00588 find_1d = p 00589 00590 end function find_1d 00591 00592 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00593 00594 !> 1D linear interpolation. 00595 !! 00596 !! Takes an array of x and y, and a single x value, and returns the 00597 !! corresponding y using linear interpolation. x_vals must be 00598 !! sorted. 00599 real(kind=dp) function interp_1d(x_vals, y_vals, x) 00600 00601 !> X value array, must be sorted. 00602 real(kind=dp), intent(in) :: x_vals(:) 00603 !> Y value array. 00604 real(kind=dp), intent(in) :: y_vals(size(x_vals)) 00605 !> Value to interpolate at. 00606 real(kind=dp), intent(in) :: x 00607 00608 integer :: n, p 00609 real(kind=dp) :: y, alpha 00610 00611 n = size(x_vals) 00612 p = find_1d(n, x_vals, x) 00613 if (p == 0) then 00614 y = y_vals(1) 00615 elseif (p == n) then 00616 y = y_vals(n) 00617 else 00618 alpha = (x - x_vals(p)) / (x_vals(p+1) - x_vals(p)) 00619 y = (1d0 - alpha) * y_vals(p) + alpha * y_vals(p+1) 00620 end if 00621 interp_1d = y 00622 00623 end function interp_1d 00624 00625 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00626 00627 !> Linear interpolation over discrete indices. 00628 !! 00629 !! Takes real values \c x_1 and \c x_n and integers \c n and \c i 00630 !! and returns the linear interpolation so that \c x_1 is returned 00631 !! when \c i = 1 and \c x_n is returned when \c i = \c n. 00632 real(kind=dp) function interp_linear_disc(x_1, x_n, n, i) 00633 00634 !> Value of \c x when \c i = 1. 00635 real(kind=dp), intent(in) :: x_1 00636 !> Value of \c x when \c i = n. 00637 real(kind=dp), intent(in) :: x_n 00638 !> Number of points to interpolate over. 00639 integer, intent(in) :: n 00640 !> Index to interpolate at. 00641 integer, intent(in) :: i 00642 00643 real(kind=dp) :: alpha 00644 00645 alpha = real(i - 1, kind=dp) / real(n - 1, kind=dp) 00646 interp_linear_disc = (1d0 - alpha) * x_1 + alpha * x_n 00647 00648 end function interp_linear_disc 00649 00650 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00651 00652 !> Convert a string to an integer. 00653 integer function string_to_integer(string) 00654 00655 !> String to convert. 00656 character(len=*), intent(in) :: string 00657 00658 integer :: val 00659 integer :: ios 00660 character(len=len(string)+300) :: error_msg 00661 00662 call assert_msg(447772570, len_trim(string) <= 20, & 00663 'error converting "' // trim(string) & 00664 // '" to integer: string too long') 00665 read(string, '(i20)', iostat=ios) val 00666 call assert_msg(895881873, ios == 0, & 00667 'error converting "' // trim(string) & 00668 // '" to integer: IOSTAT = ' // trim(integer_to_string(ios))) 00669 string_to_integer = val 00670 00671 end function string_to_integer 00672 00673 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00674 00675 !> Convert a string to a real. 00676 real(kind=dp) function string_to_real(string) 00677 00678 !> String to convert. 00679 character(len=*), intent(in) :: string 00680 00681 real(kind=dp) :: val 00682 integer :: ios 00683 character(len=len(string)+300) :: error_msg 00684 00685 call assert_msg(733728030, len_trim(string) <= 30, & 00686 'error converting "' // trim(string) // '" to real: string too long') 00687 read(string, '(f30.0)', iostat=ios) val 00688 call assert_msg(727430097, ios == 0, & 00689 'error converting "' // trim(string) & 00690 // '" to real: IOSTAT = ' // trim(integer_to_string(ios))) 00691 string_to_real = val 00692 00693 end function string_to_real 00694 00695 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00696 00697 !> Convert a string to a logical. 00698 logical function string_to_logical(string) 00699 00700 !> String to convert. 00701 character(len=*), intent(in) :: string 00702 00703 logical :: val 00704 integer :: ios 00705 character(len=len(string)+300) :: error_msg 00706 00707 val = .false. 00708 if ((trim(string) == 'yes') & 00709 .or. (trim(string) == 'y') & 00710 .or. (trim(string) == 'true') & 00711 .or. (trim(string) == 't') & 00712 .or. (trim(string) == '1')) then 00713 val = .true. 00714 elseif ((trim(string) == 'no') & 00715 .or. (trim(string) == 'n') & 00716 .or. (trim(string) == 'false') & 00717 .or. (trim(string) == 'f') & 00718 .or. (trim(string) == '0')) then 00719 val = .false. 00720 else 00721 call assert_msg(985010153, ios == 0, & 00722 'error converting "' // trim(string) & 00723 // '" to logical: IOSTAT = ' // trim(integer_to_string(ios))) 00724 end if 00725 string_to_logical = val 00726 00727 end function string_to_logical 00728 00729 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00730 00731 !> Convert an integer to a string format. 00732 character(len=PMC_UTIL_CONVERT_STRING_LEN) function integer_to_string(val) 00733 00734 !> Value to convert. 00735 integer, intent(in) :: val 00736 00737 character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val 00738 00739 ret_val = "" 00740 write(ret_val, '(i30)') val 00741 integer_to_string = adjustl(ret_val) 00742 00743 end function integer_to_string 00744 00745 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00746 00747 !> Convert a real to a string format. 00748 character(len=PMC_UTIL_CONVERT_STRING_LEN) function real_to_string(val) 00749 00750 !> Value to convert. 00751 real(kind=dp), intent(in) :: val 00752 00753 character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val 00754 00755 ret_val = "" 00756 write(ret_val, '(g30.20)') val 00757 real_to_string = adjustl(ret_val) 00758 00759 end function real_to_string 00760 00761 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00762 00763 !> Convert a logical to a string format. 00764 character(len=PMC_UTIL_CONVERT_STRING_LEN) function logical_to_string(val) 00765 00766 !> Value to convert. 00767 logical, intent(in) :: val 00768 00769 character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val 00770 00771 ret_val = "" 00772 if (val) then 00773 ret_val = "TRUE" 00774 else 00775 ret_val = "FALSE" 00776 end if 00777 logical_to_string = ret_val 00778 00779 end function logical_to_string 00780 00781 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00782 00783 !> Convert a complex to a string format. 00784 character(len=PMC_UTIL_CONVERT_STRING_LEN) function complex_to_string(val) 00785 00786 !> Value to convert. 00787 complex(kind=dc), intent(in) :: val 00788 00789 character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val 00790 00791 ret_val = "" 00792 ret_val = "(" // trim(real_to_string(real(val))) 00793 // ", " // trim(real_to_string(aimag(val))) // ")" 00794 complex_to_string = adjustl(ret_val) 00795 00796 end function complex_to_string 00797 00798 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00799 00800 !> Convert an integer to a string format of maximum length. 00801 character(len=PMC_UTIL_CONVERT_STRING_LEN) 00802 function integer_to_string_max_len(val, max_len) 00803 00804 !> Value to convert. 00805 integer, intent(in) :: val 00806 !> Maximum length of resulting string. 00807 integer, intent(in) :: max_len 00808 00809 character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val 00810 00811 ret_val = integer_to_string(val) 00812 if (len_trim(ret_val) > max_len) then 00813 ret_val = real_to_string_max_len(real(val, kind=dp), max_len) 00814 end if 00815 integer_to_string_max_len = adjustl(ret_val) 00816 00817 end function integer_to_string_max_len 00818 00819 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00820 00821 !> Convert a real to a string format of maximum length. 00822 character(len=PMC_UTIL_CONVERT_STRING_LEN) 00823 function real_to_string_max_len(val, max_len) 00824 00825 !> Value to convert. 00826 real(kind=dp), intent(in) :: val 00827 !> Maximum length of resulting string. 00828 integer, intent(in) :: max_len 00829 00830 character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val, exp_str, frac_str 00831 integer :: exp_val, exp_len, frac_len, use_frac_len, min_frac_len, i 00832 real(kind=dp) :: frac_val 00833 00834 ret_val = "" 00835 if (val == 0d0) then 00836 if (max_len >= 3) then 00837 ret_val = "0e0" 00838 else 00839 do i = 1,max_len 00840 ret_val(i:i) = "*" 00841 end do 00842 end if 00843 real_to_string_max_len = adjustl(ret_val) 00844 return 00845 end if 00846 00847 exp_val = floor(log10(abs(val))) 00848 frac_val = val / 10d0**exp_val 00849 exp_str = integer_to_string(exp_val) 00850 frac_str = real_to_string(frac_val) 00851 00852 exp_len = len_trim(exp_str) 00853 frac_len = len_trim(frac_str) 00854 use_frac_len = max_len - 1 - exp_len 00855 if (use_frac_len > frac_len) then 00856 use_frac_len = frac_len 00857 end if 00858 if (val < 0d0) then 00859 min_frac_len = 2 00860 else 00861 min_frac_len = 1 00862 end if 00863 if (use_frac_len < min_frac_len) then 00864 do i = 1,max_len 00865 ret_val(i:i) = "*" 00866 end do 00867 else 00868 ret_val = frac_str(1:use_frac_len) // "e" // trim(exp_str) 00869 end if 00870 real_to_string_max_len = adjustl(ret_val) 00871 00872 end function real_to_string_max_len 00873 00874 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00875 00876 !> Convert a time to a string format of maximum length. 00877 character(len=PMC_UTIL_CONVERT_STRING_LEN) 00878 function time_to_string_max_len(time, max_len) 00879 00880 !> Time to convert (s). 00881 real(kind=dp), intent(in) :: time 00882 !> Maximum length of resulting string. 00883 integer, intent(in) :: max_len 00884 00885 integer, dimension(4), parameter :: scale = (/ 1, 60, 60, 24 /) 00886 character, dimension(4), parameter :: unit = (/ "s", "m", "h", "d" /) 00887 00888 character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val 00889 integer :: i 00890 logical :: len_ok 00891 real(kind=dp) :: scaled_time 00892 00893 scaled_time = time 00894 len_ok = .false. 00895 do i = 1,4 00896 scaled_time = scaled_time / real(scale(i), kind=dp) 00897 ret_val = trim(integer_to_string(nint(scaled_time))) // unit(i) 00898 if (len_trim(ret_val) <= max_len) then 00899 len_ok = .true. 00900 exit 00901 end if 00902 end do 00903 if (.not. len_ok) then 00904 ret_val = "" 00905 do i = 1,max_len 00906 ret_val(i:i) = "*" 00907 end do 00908 end if 00909 time_to_string_max_len = adjustl(ret_val) 00910 00911 end function time_to_string_max_len 00912 00913 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00914 00915 !> Convert a real-valued vector into an integer-valued vector. 00916 !! 00917 !! Use n_samp samples. Returns discrete vector whose relative entry 00918 !! sizes are \f$ \ell_1 \f$ closest to the continuous vector. 00919 subroutine vec_cts_to_disc(n, vec_cts, n_samp, vec_disc) 00920 00921 !> Number of entries in vectors. 00922 integer, intent(in) :: n 00923 !> Continuous vector. 00924 real(kind=dp), intent(in) :: vec_cts(n) 00925 !> Number of discrete samples to use. 00926 integer, intent(in) :: n_samp 00927 !> Discretized vector. 00928 integer, intent(out) :: vec_disc(n) 00929 00930 integer :: k(1) 00931 real(kind=dp) :: vec_tot 00932 00933 vec_tot = sum(vec_cts) 00934 00935 ! assign a best guess for each bin independently 00936 vec_disc = nint(vec_cts / vec_tot * real(n_samp, kind=dp)) 00937 00938 ! if we have too few particles then add more 00939 do while (sum(vec_disc) < n_samp) 00940 k = minloc(abs(real(vec_disc + 1, kind=dp) - vec_cts) 00941 - abs(real(vec_disc, kind=dp) - vec_cts)) 00942 vec_disc(k) = vec_disc(k) + 1 00943 end do 00944 00945 ! if we have too many particles then remove some 00946 do while (sum(vec_disc) > n_samp) 00947 k = minloc(abs(real(vec_disc - 1, kind=dp) - vec_cts) 00948 - abs(real(vec_disc, kind=dp) - vec_cts)) 00949 vec_disc(k) = vec_disc(k) - 1 00950 end do 00951 00952 call assert_msg(323412496, sum(vec_disc) == n_samp, & 00953 'generated incorrect number of samples') 00954 00955 end subroutine vec_cts_to_disc 00956 00957 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00958 00959 !> Computes the average of an array of integer numbers. 00960 subroutine average_integer(int_vec, int_avg) 00961 00962 !> Array of integer numbers. 00963 integer, intent(in) :: int_vec(:) 00964 !> Average of int_vec. 00965 integer, intent(out) :: int_avg 00966 00967 int_avg = sum(int_vec) / size(int_vec) 00968 00969 end subroutine average_integer 00970 00971 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00972 00973 !> Computes the average of an array of real numbers. 00974 subroutine average_real(real_vec, real_avg) 00975 00976 !> Array of real numbers. 00977 real(kind=dp), intent(in) :: real_vec(:) 00978 !> Average of real_vec. 00979 real(kind=dp), intent(out) :: real_avg 00980 00981 real_avg = sum(real_vec) / real(size(real_vec), kind=dp) 00982 00983 end subroutine average_real 00984 00985 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00986 00987 !> Reallocate the given array to ensure it is of the given size. 00988 subroutine ensure_real_array_size(x, n) 00989 00990 !> Array of real numbers. 00991 real(kind=dp), intent(inout), allocatable :: x(:) 00992 !> Desired size of array. 00993 integer, intent(in) :: n 00994 00995 if (size(x) /= n) then 00996 deallocate(x) 00997 allocate(x(n)) 00998 end if 00999 01000 end subroutine ensure_real_array_size 01001 01002 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01003 01004 !> Strip the extension to find the basename. 01005 !! 01006 !! The filename is assumed to be of the form basename.extension, 01007 !! where extension contains no periods. If there is no period in 01008 !! filename then basename is the whole filename. 01009 subroutine get_basename(filename, basename) 01010 01011 !> Full filename. 01012 character(len=*), intent(in) :: filename 01013 !> Basename. 01014 character(len=*), intent(out) :: basename 01015 01016 integer :: i 01017 logical :: found_period 01018 01019 basename = filename 01020 i = len_trim(basename) 01021 found_period = .false. 01022 do while ((i > 0) .and. (.not. found_period)) 01023 ! Fortran .and. does not short-circuit, so we can't do the 01024 ! obvious do while ((i > 0) .and. (basename(i:i) /= ".")), 01025 ! instead we have to use this hack with the found_period 01026 ! logical variable. 01027 if (basename(i:i) == ".") then 01028 found_period = .true. 01029 else 01030 i = i - 1 01031 end if 01032 end do 01033 if (i > 0) then 01034 basename(i:) = "" 01035 end if 01036 01037 end subroutine get_basename 01038 01039 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01040 01041 !> Current date and time in ISO 8601 format. 01042 subroutine iso8601_date_and_time(date_time) 01043 01044 !> Result string. 01045 character(len=*), intent(out) :: date_time 01046 01047 character(len=10) :: date, time, zone 01048 01049 call assert_msg(893219839, len(date_time) >= 29, & 01050 "date_time string must have length at least 29") 01051 call date_and_time(date, time, zone) 01052 date_time = "" 01053 write(date_time, '(14a)') date(1:4), "-", date(5:6), "-", & 01054 date(7:8), "T", time(1:2), ":", time(3:4), ":", & 01055 time(5:10), zone(1:3), ":", zone(4:5) 01056 01057 end subroutine iso8601_date_and_time 01058 01059 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01060 01061 !> Convert degrees to radians. 01062 real(kind=dp) function deg2rad(deg) 01063 01064 !> Input degrees. 01065 real(kind=dp), intent(in) :: deg 01066 01067 deg2rad = deg / 180d0 * const%pi 01068 01069 end function deg2rad 01070 01071 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01072 01073 !> Convert radians to degrees. 01074 real(kind=dp) function rad2deg(rad) 01075 01076 !> Input radians. 01077 real(kind=dp), intent(in) :: rad 01078 01079 rad2deg = rad / const%pi * 180d0 01080 01081 end function rad2deg 01082 01083 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01084 01085 #ifdef DEFINE_LOCAL_COMMAND_ARGUMENT 01086 !> Hack for compilers that don't support the Fortran2003 standard 01087 !> command_argument_count() function. 01088 integer function command_argument_count() 01089 01090 command_argument_count = iargc() 01091 01092 end function command_argument_count 01093 #endif 01094 01095 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01096 01097 #ifdef DEFINE_LOCAL_COMMAND_ARGUMENT 01098 !> Hack for compilers that don't support the Fortran2003 standard 01099 !> get_command_argument() subroutine. 01100 subroutine get_command_argument(i, arg) 01101 integer, intent(in) :: i 01102 character(len=*), intent(out) :: arg 01103 01104 call getarg(i, arg) 01105 01106 end subroutine get_command_argument 01107 #endif 01108 01109 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01110 01111 !> Sort the given data array and return the permutation defining the sort. 01112 !! 01113 !! If \c data is the original data array, \c new_data is the sorted 01114 !! value of data, and \c perm is the returned permutation, then 01115 !! <tt>new_data(i) = data(perm(i))</tt>. 01116 subroutine integer_sort(data, perm) 01117 01118 !> Data array to sort, sorted on exit. 01119 integer, intent(inout) :: data(:) 01120 !> Permutation defining the sort: <tt>new_data(i) = data(perm(i))</tt>. 01121 integer, intent(out) :: perm(size(data)) 01122 01123 #ifdef PMC_USE_C_SORT 01124 integer(kind=c_int) :: n_c 01125 integer(kind=c_int), target :: data_c(size(data)) 01126 integer(kind=c_int), target :: perm_c(size(data)) 01127 type(c_ptr) :: data_ptr, perm_ptr 01128 01129 #ifndef DOXYGEN_SKIP_DOC 01130 interface 01131 subroutine integer_sort_c(n_c, data_ptr, perm_ptr) bind(c) 01132 use iso_c_binding 01133 integer(kind=c_int), value :: n_c 01134 type(c_ptr), value :: data_ptr, perm_ptr 01135 end subroutine integer_sort_c 01136 end interface 01137 #endif 01138 01139 data_c = int(data, kind=c_int) 01140 perm_c = 0_c_int 01141 n_c = int(size(data), kind=c_int) 01142 data_ptr = c_loc(data_c) 01143 perm_ptr = c_loc(perm_c) 01144 call integer_sort_c(n_c, data_ptr, perm_ptr) 01145 data = int(data_c) 01146 perm = int(perm_c) 01147 #endif 01148 01149 end subroutine integer_sort 01150 01151 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01152 01153 !> Load a real array from a text file. 01154 subroutine loadtxt(filename, data) 01155 01156 !> Filename to read from. 01157 character(len=*), intent(in) :: filename 01158 !> Array to store data in. 01159 real(kind=dp), intent(inout), allocatable :: data(:,:) 01160 01161 integer :: unit, row, col 01162 logical :: eol, eof 01163 character(len=1000) :: word 01164 01165 deallocate(data) 01166 allocate(data(1,0)) 01167 call open_file_read(filename, unit) 01168 01169 eof = .false. 01170 row = 1 01171 col = 1 01172 do while (.not. eof) 01173 call read_word_raw(unit, word, eol, eof) 01174 if (len_trim(word) > 0) then 01175 if (row == 1) then 01176 if (col > size(data, 2)) then 01177 call reallocate_real_array2d(data, 1, 2 * col) 01178 end if 01179 else 01180 if (col > size(data, 2)) then 01181 call assert_msg(516120334, col <= size(data, 2), & 01182 trim(filename) // ": line " & 01183 // trim(integer_to_string(row)) // " too long") 01184 end if 01185 end if 01186 if (row > size(data, 1)) then 01187 call reallocate_real_array2d(data, 2 * row, size(data, 2)) 01188 end if 01189 data(row, col) = string_to_real(word) 01190 col = col + 1 01191 end if 01192 if (eol .or. eof) then 01193 if (row == 1) then 01194 call reallocate_real_array2d(data, 1, col - 1) 01195 else 01196 call assert_msg(266924956, & 01197 (col == 1) .or. (col == size(data, 2) + 1), & 01198 trim(filename) // ": line " & 01199 // trim(integer_to_string(row)) // " too short") 01200 end if 01201 end if 01202 if (eol) then 01203 row = row + 1 01204 col = 1 01205 end if 01206 end do 01207 01208 if (col == 1) then 01209 call reallocate_real_array2d(data, row - 1, size(data, 2)) 01210 end if 01211 call close_file(unit) 01212 01213 end subroutine loadtxt 01214 01215 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01216 01217 !> Reallocate a 2D real array to the given size, preserving the contents. 01218 subroutine reallocate_real_array2d(data, rows, cols) 01219 01220 !> Array to reallocate. 01221 real(kind=dp), allocatable, intent(inout) :: data(:,:) 01222 !> New leading dimension. 01223 integer, intent(in) :: rows 01224 !> New trailing dimension. 01225 integer, intent(in) :: cols 01226 01227 real(kind=dp) :: tmp_data(rows, cols) 01228 integer :: data_rows, data_cols 01229 01230 data_rows = min(rows, size(data, 1)) 01231 data_cols = min(cols, size(data, 2)) 01232 tmp_data(1:data_rows, 1:data_cols) = data(1:data_rows, 1:data_cols) 01233 deallocate(data) 01234 allocate(data(rows, cols)) 01235 data(1:data_rows, 1:data_cols) = tmp_data(1:data_rows, 1:data_cols) 01236 01237 end subroutine reallocate_real_array2d 01238 01239 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01240 01241 !> Read a single character from a file, signaling if we have hit 01242 !> end-of-line (EOL) or end-of-file (EOF). If EOL or EOF are true 01243 !> then the character value should be ignored. 01244 !! 01245 !! Testing with gfortran 4.5.3 and different length files shows: 01246 !! 01247 !! Empty file (total file length 0): 01248 !! * eol = .false., eof = .true. 01249 !! * subsequent calls error 01250 !! 01251 !! File containing a single 'A' character (total file length 1): 01252 !! * char = 'A', eol = .false., eof = .false. 01253 !! * eol = .false., eof = .true. 01254 !! * subsequent calls error 01255 !! 01256 !! File containing a single newline '\n' (total file length 1): 01257 !! * eol = .true., eof = .false. 01258 !! * eol = .false., eof = .true. 01259 !! * subsequent calls error 01260 !! 01261 !! File containing a character and newline 'A\n' (total file length 2): 01262 !! * char = 'A', eol = .false., eof = .false. 01263 !! * eol = .true., eof = .false. 01264 !! * eol = .false., eof = .true. 01265 !! * subsequent calls error 01266 subroutine read_char_raw(unit, char, eol, eof) 01267 01268 !> Unit number to read from. 01269 integer, intent(in) :: unit 01270 !> Character read. 01271 character, intent(out) :: char 01272 !> True if at EOL (end of line). 01273 logical, intent(out) :: eol 01274 !> True if at EOF (end of file). 01275 logical, intent(out) :: eof 01276 01277 integer :: ios, n_read 01278 character(len=1) :: read_char 01279 01280 eol = .false. 01281 eof = .false. 01282 char = " " ! shut up uninitialized variable warnings 01283 read_char = "" ! needed for pgf95 for reading blank lines 01284 read(unit=unit, fmt='(a)', advance='no', end=100, eor=110, & 01285 iostat=ios) read_char 01286 if (ios /= 0) then 01287 write(0,*) 'ERROR: reading file: IOSTAT = ', ios 01288 stop 2 01289 end if 01290 ! only reach here if we didn't hit end-of-record (end-of-line) in 01291 ! the above read 01292 char = read_char 01293 goto 120 01294 01295 100 eof = .true. ! goto here if end-of-file was encountered immediately 01296 goto 120 01297 01298 110 eol = .true. ! goto here if end-of-record, meaning end-of-line 01299 01300 120 return 01301 01302 end subroutine read_char_raw 01303 01304 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01305 01306 !> Read a white-space delimited word from a file, signaling if we 01307 !> have EOL or EOF. If EOL or EOF are true then the word will still 01308 !> be meaningful data. If there was no data to be read then 01309 !> len(word) will be 0. 01310 subroutine read_word_raw(unit, word, eol, eof) 01311 01312 !> Unit number to read from. 01313 integer, intent(in) :: unit 01314 !> Word read. 01315 character(len=*), intent(out) :: word 01316 !> True if at EOL (end of line). 01317 logical, intent(out) :: eol 01318 !> True if at EOF (end of file). 01319 logical, intent(out) :: eof 01320 01321 integer :: i 01322 character :: char 01323 01324 word = "" 01325 01326 ! skip over spaces 01327 call read_char_raw(unit, char, eol, eof) 01328 do while (((ichar(char) == 9) .or. (ichar(char) == 32)) & 01329 .and. (.not. eol) .and. (.not. eof)) 01330 call read_char_raw(unit, char, eol, eof) 01331 end do 01332 if (eol .or. eof) return 01333 01334 ! char is now the first word character 01335 i = 1 01336 word(i:i) = char 01337 call read_char_raw(unit, char, eol, eof) 01338 do while ((ichar(char) /= 9) .and. (ichar(char) /= 32) & 01339 .and. (.not. eol) .and. (.not. eof) .and. (i < len(word))) 01340 i = i + 1 01341 word(i:i) = char 01342 if (i < len(word)) then 01343 call read_char_raw(unit, char, eol, eof) 01344 end if 01345 end do 01346 01347 end subroutine read_word_raw 01348 01349 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01350 01351 end module pmc_util