PartMC  2.3.0
util.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2015 Nicole Riemer and Matthew West
2 ! Licensed under the GNU General Public License version 2 or (at your
3 ! option) any later version. See the file COPYING for details.
4 
5 !> \file
6 !> The pmc_util module.
7 
8 !> Common utility subroutines.
9 module pmc_util
10 
11  use pmc_constants
12 #ifdef PMC_USE_MPI
13  use mpi
14 #endif
15 #ifdef PMC_USE_C_SORT
16  use iso_c_binding
17 #endif
18 
19  !> Maximum number of IO units usable simultaneously.
20  integer, parameter :: max_units = 200
21  !> Minimum unit number to allocate.
22  integer, parameter :: unit_offset = 10
23  !> Table of unit numbers storing allocation status.
24  logical, save :: unit_used(max_units) = .false.
25 
26  !> Length of string for converting numbers.
27  integer, parameter :: PMC_UTIL_CONVERT_STRING_LEN = 100
28  !> Maximum length of filenames.
29  integer, parameter :: PMC_MAX_FILENAME_LEN = 300
30 
31 contains
32 
33 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
34 
35  !> Prints a warning message.
36  subroutine warn_msg(code, warning_msg, already_warned)
37 
38  !> Status code to use.
39  integer, intent(in) :: code
40  !> Message to display.
41  character(len=*), intent(in) :: warning_msg
42  !> Flag to control warning only once (should be a save variable).
43  logical, intent(inout), optional :: already_warned
44 
45  if (present(already_warned)) then
46  if (already_warned) return
47  ! set already_warned so next time we will immediately return
48  already_warned = .true.
49  end if
50  write(0,'(a)') 'WARNING (PartMC-' // trim(integer_to_string(code)) &
51  // '): ' // trim(warning_msg)
52 
53  end subroutine warn_msg
54 
55 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
56 
57  !> Prints a warning message if condition_ok is false.
58  subroutine warn_assert_msg(code, condition_ok, warning_msg)
59 
60  !> Status code to use.
61  integer, intent(in) :: code
62  !> Whether the assertion is ok.
63  logical, intent(in) :: condition_ok
64  !> Message to display.
65  character(len=*), intent(in) :: warning_msg
66 
67  if (.not. condition_ok) then
68  call warn_msg(code, warning_msg)
69  end if
70 
71  end subroutine warn_assert_msg
72 
73 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
74 
75  !> Errors unless condition_ok is true.
76  subroutine assert_msg(code, condition_ok, error_msg)
77 
78  !> Status code to use if assertion fails.
79  integer, intent(in) :: code
80  !> Whether the assertion is ok.
81  logical, intent(in) :: condition_ok
82  !> Msg if assertion fails.
83  character(len=*), intent(in) :: error_msg
84 
85  integer :: ierr
86 
87  if (.not. condition_ok) then
88  write(0,'(a)') 'ERROR (PartMC-' // trim(integer_to_string(code)) &
89  // '): ' // trim(error_msg)
90 #ifdef PMC_USE_MPI
91  call mpi_abort(mpi_comm_world, code, ierr)
92 #else
93  stop 3
94 #endif
95  end if
96 
97  end subroutine assert_msg
98 
99 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
100 
101  !> Errors unless condition_ok is true.
102  subroutine assert(code, condition_ok)
103 
104  !> Status code to use if assertion fails.
105  integer, intent(in) :: code
106  !> Whether the assertion is ok.
107  logical, intent(in) :: condition_ok
108 
109  if (.not. condition_ok) then
110  ! much faster with gfortran 4.4.5 to do this extra test
111  ! FIXME: is it still faster now that assert_msg doesn't
112  ! unconditionally construct a code_str?
113  call assert_msg(code, condition_ok, 'assertion failed')
114  end if
115 
116  end subroutine assert
117 
118 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
119 
120  !> Error immediately.
121  subroutine die(code)
122 
123  !> Status code to use if assertion fails.
124  integer, intent(in) :: code
125 
126  call assert(code, .false.)
127 
128  end subroutine die
129 
130 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
131 
132  !> Error immediately.
133  subroutine die_msg(code, error_msg)
134 
135  !> Status code to use if assertion fails.
136  integer, intent(in) :: code
137  !> Msg if assertion fails.
138  character(len=*), intent(in) :: error_msg
139 
140  call assert_msg(code, .false., error_msg)
141 
142  end subroutine die_msg
143 
144 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
145 
146  !> Returns an available unit number. This should be freed by free_unit().
147  integer function get_unit()
148 
149  integer i
150  logical found_unit
151 
152  found_unit = .false.
153  do i = 1,max_units
154  if (.not. unit_used(i)) then
155  found_unit = .true.
156  exit
157  end if
158  end do
159  if (.not. found_unit) then
160  call die_msg(690355443, &
161  'no more units available - need to free_unit()')
162  end if
163  unit_used(i) = .true.
164  get_unit = i + unit_offset
165 
166  end function get_unit
167 
168 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
169 
170  !> Frees a unit number returned by get_unit().
171  subroutine free_unit(unit)
172 
173  integer, intent(in) :: unit
174 
175  unit_used(unit - unit_offset) = .false.
176 
177  end subroutine free_unit
178 
179 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
180 
181  !> Open a file for reading with an automatically assigned unit and
182  !> test that it succeeds. The file should be closed with
183  !> close_file().
184  subroutine open_file_read(filename, unit)
185 
186  !> Filename to open.
187  character(len=*), intent(in) :: filename
188  !> Unit assigned to file.
189  integer, intent(out) :: unit
190 
191  integer :: ios
192 
193  unit = get_unit()
194  open(unit=unit, file=filename, status='old', action='read', iostat=ios)
195  call assert_msg(544344918, ios == 0, 'unable to open file ' &
196  // trim(filename) // ' for reading: ' // integer_to_string(ios))
197 
198  end subroutine open_file_read
199 
200 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
201 
202  !> Open a file for writing with an automatically assigned unit and
203  !> test that it succeeds. The file should be closed with
204  !> close_file().
205  subroutine open_file_write(filename, unit)
206 
207  !> Filename to open.
208  character(len=*), intent(in) :: filename
209  !> Unit assigned to file.
210  integer, intent(out) :: unit
211 
212  integer :: ios
213 
214  unit = get_unit()
215  open(unit=unit, file=filename, status='replace', action='write', &
216  iostat=ios)
217  call assert_msg(609624199, ios == 0, 'unable to open file ' &
218  // trim(filename) // ' for writing: ' // integer_to_string(ios))
219 
220  end subroutine open_file_write
221 
222 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
223 
224  !> Close a file and de-assign the unit.
225  subroutine close_file(unit)
226 
227  !> Unit to close.
228  integer, intent(in) :: unit
229 
230  close(unit)
231  call free_unit(unit)
232 
233  end subroutine close_file
234 
235 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
236 
237  !> Convert volume (m^3) to radius (m).
238  real(kind=dp) elemental function vol2rad(v)
239 
240  !> Volume (m^3).
241  real(kind=dp), intent(in) :: v
242 
243  vol2rad = (v / (4d0 / 3d0 * const%pi))**(1d0/3d0)
244 
245  end function vol2rad
246 
247 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
248 
249  !> Convert volume (m^3) to diameter (m).
250  real(kind=dp) elemental function vol2diam(v)
251 
252  !> Volume (m^3).
253  real(kind=dp), intent(in) :: v
254 
255  vol2diam = 2d0 * (v / (4d0 / 3d0 * const%pi))**(1d0/3d0)
256 
257  end function vol2diam
258 
259 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
260 
261  !> Convert radius (m) to diameter (m).
262  real(kind=dp) elemental function rad2diam(r)
263 
264  !> Radius (m).
265  real(kind=dp), intent(in) :: r
266 
267  rad2diam = 2d0 * r
268 
269  end function rad2diam
270 
271 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
272 
273  !> Convert radius (m) to volume (m^3).
274  real(kind=dp) elemental function rad2vol(r)
275 
276  !> Radius (m).
277  real(kind=dp), intent(in) :: r
278 
279  rad2vol = 4d0 / 3d0 * const%pi * r**3d0
280 
281  end function rad2vol
282 
283 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
284 
285  !> Convert diameter (m) to radius (m).
286  real(kind=dp) elemental function diam2rad(d)
287 
288  !> Diameter (m).
289  real(kind=dp), intent(in) :: d
290 
291  diam2rad = d / 2d0
292 
293  end function diam2rad
294 
295 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
296 
297  !> Convert diameter (m) to volume (m^3).
298  real(kind=dp) elemental function diam2vol(d)
299 
300  !> Diameter (m).
301  real(kind=dp), intent(in) :: d
302 
303  diam2vol = 4d0 / 3d0 * const%pi * (d / 2d0)**3d0
304 
305  end function diam2vol
306 
307 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
308 
309  !> Tests whether two real numbers are almost equal using only a
310  !> relative tolerance.
311  logical function almost_equal(d1, d2)
312 
313  !> First number to compare.
314  real(kind=dp), intent(in) :: d1
315  !> Second number to compare.
316  real(kind=dp), intent(in) :: d2
317 
318  !> Relative tolerance.
319  real(kind=dp), parameter :: eps = 1d-8
320 
321  ! handle the 0.0 case
322  if (d1 .eq. d2) then
323  almost_equal = .true.
324  else
325  if (abs(d1 - d2) / (abs(d1) + abs(d2)) .lt. eps) then
326  almost_equal = .true.
327  else
328  almost_equal = .false.
329  end if
330  end if
331 
332  end function almost_equal
333 
334 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
335 
336  !> Tests whether two real numbers are almost equal using an
337  !> absolute and relative tolerance.
338  logical function almost_equal_abs(d1, d2, abs_tol)
339 
340  !> First number to compare.
341  real(kind=dp), intent(in) :: d1
342  !> Second number to compare.
343  real(kind=dp), intent(in) :: d2
344  !> Tolerance for when d1 equals d2.
345  real(kind=dp), intent(in) :: abs_tol
346 
347  !> Relative tolerance.
348  real(kind=dp), parameter :: eps = 1d-8
349 
350  ! handle the 0.0 case
351  if (d1 .eq. d2) then
352  almost_equal_abs = .true.
353  else
354  if ((abs(d1 - d2) .lt. abs_tol) .or. &
355  (abs(d1 - d2) / (abs(d1) + abs(d2)) .lt. eps)) then
356  almost_equal_abs = .true.
357  else
358  almost_equal_abs = .false.
359  end if
360  end if
361 
362  end function almost_equal_abs
363 
364 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
365 
366  !> Check that the first time interval is close to an integer
367  !> multiple of the second, and warn if it is not.
368  subroutine check_time_multiple(first_name, first_time, &
369  second_name, second_time)
370 
371  !> Name of the first time variable.
372  character(len=*), intent(in) :: first_name
373  !> First time variable (s).
374  real(kind=dp), intent(in) :: first_time
375  !> Name of the second time variable.
376  character(len=*), intent(in) :: second_name
377  !> Second time variable (s).
378  real(kind=dp), intent(in) :: second_time
379 
380  real(kind=dp) :: ratio
381 
382  ratio = first_time / second_time
383  if (abs(ratio - aint(ratio)) > 1d-6) then
384  call warn_msg(952299377, trim(first_name) &
385  // " is not an integer multiple of " // trim(second_name))
386  end if
387 
388  end subroutine check_time_multiple
389 
390 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
391 
392  !> Computes whether an event is scheduled to take place.
393  !!
394  !! The events should occur ideally at times 0, interval, 2*interval,
395  !! etc. The events are guaranteed to occur at least interval * (1 -
396  !! tolerance) apart, and if at least interval time has passed then
397  !! the next call is guaranteed to do the event. Otherwise the
398  !! timestep is used to guess whether to do the event.
399  subroutine check_event(time, timestep, interval, last_time, &
400  do_event)
401 
402  !> Current time.
403  real(kind=dp), intent(in) :: time
404  !> Estimate of the time to the next call.
405  real(kind=dp), intent(in) :: timestep
406  !> How often the event should be done.
407  real(kind=dp), intent(in) :: interval
408  !> When the event was last done.
409  real(kind=dp), intent(inout) :: last_time
410  !> Whether the event should be done.
411  logical, intent(out) :: do_event
412 
413  !> Fuzz for event occurance.
414  real(kind=dp), parameter :: tolerance = 1d-6
415 
416  real(kind=dp) closest_interval_time
417 
418  ! if we are at time 0 then do the event unconditionally
419  if (time .eq. 0d0) then
420  do_event = .true.
421  else
422  ! if we are too close to the last time then don't do it
423  if ((time - last_time) .lt. interval * (1d0 - tolerance)) then
424  do_event = .false.
425  else
426  ! if it's been too long since the last time then do it
427  if ((time - last_time) .ge. interval) then
428  do_event = .true.
429  else
430  ! gray area -- if we are closer than we will be next
431  ! time then do it
432  closest_interval_time = anint(time / interval) * interval
433  if (abs(time - closest_interval_time) &
434  .lt. abs(time + timestep - closest_interval_time)) &
435  then
436  do_event = .true.
437  else
438  do_event = .false.
439  end if
440  end if
441  end if
442  end if
443 
444  if (do_event) then
445  last_time = time
446  end if
447 
448  end subroutine check_event
449 
450 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
451 
452  !> Makes a linearly spaced array from min to max.
453  subroutine linspace(min_x, max_x, x)
454 
455  !> Minimum array value.
456  real(kind=dp), intent(in) :: min_x
457  !> Maximum array value.
458  real(kind=dp), intent(in) :: max_x
459  !> Array.
460  real(kind=dp), intent(out) :: x(:)
461 
462  integer :: i, n
463  real(kind=dp) :: a
464 
465  n = size(x)
466  do i = 2, (n - 1)
467  a = real(i - 1, kind=dp) / real(n - 1, kind=dp)
468  x(i) = (1d0 - a) * min_x + a * max_x
469  end do
470  if (n > 0) then
471  ! make sure these values are exact
472  x(1) = min_x
473  x(n) = max_x
474  end if
475 
476  end subroutine linspace
477 
478 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
479 
480  !> Makes a logarithmically spaced array of length n from min to max.
481  subroutine logspace(min_x, max_x, x)
482 
483  !> Minimum array value.
484  real(kind=dp), intent(in) :: min_x
485  !> Maximum array value.
486  real(kind=dp), intent(in) :: max_x
487  !> Array.
488  real(kind=dp), intent(out) :: x(:)
489 
490  integer :: n
491  real(kind=dp) :: log_x(size(x))
492 
493  n = size(x)
494  if (n == 0) return
495  call assert(548290438, min_x > 0d0)
496  call assert(805259035, max_x > 0d0)
497  call linspace(log(min_x), log(max_x), log_x)
498  x = exp(log_x)
499  ! make sure these values are exact
500  x(1) = min_x
501  x(n) = max_x
502 
503  end subroutine logspace
504 
505 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
506 
507  !> Find the position of a real number in a 1D linear array.
508  !!
509  !! If xa is the array allocated by linspace(min_x, max_x, xa) then i
510  !! = linspace_find(min_x, max_x, n, x) returns the index i
511  !! satisfying xa(i) <= x < xa(i+1) for min_x <= x < max_x. If x ==
512  !! max_x then i = n - 1. If x > max_x then i = n. If x < min_x then
513  !! i = 0. Thus 0 <= i <= n. Here n is the length of xa.
514  !!
515  !! This is equivalent to using find_1d() but much faster if the
516  !! array is linear.
517  integer function linspace_find(min_x, max_x, n, x)
518 
519  !> Minimum array value.
520  real(kind=dp), intent(in) :: min_x
521  !> Maximum array value.
522  real(kind=dp), intent(in) :: max_x
523  !> Number of entries.
524  integer, intent(in) :: n
525  !> Value.
526  real(kind=dp), intent(in) :: x
527 
528  if (x == max_x) then
529  linspace_find = n - 1
530  return
531  end if
532  linspace_find = floor((x - min_x) / (max_x - min_x) &
533  * real(n - 1, kind=dp)) + 1
534  linspace_find = min(linspace_find, n)
535  linspace_find = max(linspace_find, 0)
536 
537  end function linspace_find
538 
539 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
540 
541  !> Find the position of a real number in a 1D logarithmic array.
542  !!
543  !! If xa is the array allocated by logspace(min_x, max_x, xa) then i
544  !! = logspace_find(min_x, max_x, n, x) returns the index i
545  !! satisfying xa(i) <= x < xa(i+1) for min_x <= x < max_x. If x >=
546  !! max_x then i = n. If x < min_x then i = 0. Thus 0 <= i <=
547  !! n. Here n is the length of xa.
548  !!
549  !! This is equivalent to using find_1d() but much faster if the
550  !! array is logarithmic.
551  integer function logspace_find(min_x, max_x, n, x)
552 
553  !> Minimum array value.
554  real(kind=dp), intent(in) :: min_x
555  !> Maximum array value.
556  real(kind=dp), intent(in) :: max_x
557  !> Number of entries.
558  integer, intent(in) :: n
559  !> Value.
560  real(kind=dp), intent(in) :: x
561 
562  logspace_find = linspace_find(log(min_x), log(max_x), n, log(x))
563 
564  end function logspace_find
565 
566 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
567 
568  !> Find the position of a real number in an arbitrary 1D array.
569  !!
570  !! Takes an array of x_vals, and a single x value, and returns the
571  !! position p such that x_vals(p) <= x < x_vals(p+1). If p == 0 then
572  !! x < x_vals(1) and if p == n then x_vals(n) <= x. x_vals must be
573  !! sorted in increasing order.
574  !!
575  !! If the array is known to be linearly or logarithmically spaced
576  !! then linspace_find() or logspace_find() will be much faster.
577  integer function find_1d(n, x_vals, x)
578 
579  !> Number of values.
580  integer, intent(in) :: n
581  !> X value array, must be sorted.
582  real(kind=dp), intent(in) :: x_vals(n)
583  !> Value to interpolate at.
584  real(kind=dp), intent(in) :: x
585 
586  integer p
587 
588  if (n == 0) then
589  find_1d = 0
590  return
591  end if
592  p = 1
593  do while (x >= x_vals(p))
594  p = p + 1
595  if (p > n) then
596  exit
597  end if
598  end do
599  p = p - 1
600  find_1d = p
601 
602  end function find_1d
603 
604 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
605 
606  !> 1D linear interpolation.
607  !!
608  !! Takes an array of x and y, and a single x value, and returns the
609  !! corresponding y using linear interpolation. x_vals must be
610  !! sorted.
611  real(kind=dp) function interp_1d(x_vals, y_vals, x)
612 
613  !> X value array, must be sorted.
614  real(kind=dp), intent(in) :: x_vals(:)
615  !> Y value array.
616  real(kind=dp), intent(in) :: y_vals(size(x_vals))
617  !> Value to interpolate at.
618  real(kind=dp), intent(in) :: x
619 
620  integer :: n, p
621  real(kind=dp) :: y, alpha
622 
623  n = size(x_vals)
624  p = find_1d(n, x_vals, x)
625  if (p == 0) then
626  y = y_vals(1)
627  elseif (p == n) then
628  y = y_vals(n)
629  else
630  alpha = (x - x_vals(p)) / (x_vals(p+1) - x_vals(p))
631  y = (1d0 - alpha) * y_vals(p) + alpha * y_vals(p+1)
632  end if
633  interp_1d = y
634 
635  end function interp_1d
636 
637 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
638 
639  !> Linear interpolation over discrete indices.
640  !!
641  !! Takes real values \c x_1 and \c x_n and integers \c n and \c i
642  !! and returns the linear interpolation so that \c x_1 is returned
643  !! when \c i = 1 and \c x_n is returned when \c i = \c n.
644  real(kind=dp) function interp_linear_disc(x_1, x_n, n, i)
645 
646  !> Value of \c x when \c i = 1.
647  real(kind=dp), intent(in) :: x_1
648  !> Value of \c x when \c i = n.
649  real(kind=dp), intent(in) :: x_n
650  !> Number of points to interpolate over.
651  integer, intent(in) :: n
652  !> Index to interpolate at.
653  integer, intent(in) :: i
654 
655  real(kind=dp) :: alpha
656 
657  alpha = real(i - 1, kind=dp) / real(n - 1, kind=dp)
658  interp_linear_disc = (1d0 - alpha) * x_1 + alpha * x_n
659 
660  end function interp_linear_disc
661 
662 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
663 
664  !> Convert a string to an integer.
665  integer function string_to_integer(string)
666 
667  !> String to convert.
668  character(len=*), intent(in) :: string
669 
670  integer :: val
671  integer :: ios
672  character(len=len(string)+300) :: error_msg
673 
674  call assert_msg(447772570, len_trim(string) <= 20, &
675  'error converting "' // trim(string) &
676  // '" to integer: string too long')
677  read(string, '(i20)', iostat=ios) val
678  call assert_msg(895881873, ios == 0, &
679  'error converting "' // trim(string) &
680  // '" to integer: IOSTAT = ' // trim(integer_to_string(ios)))
681  string_to_integer = val
682 
683  end function string_to_integer
684 
685 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
686 
687  !> Convert a string to a real.
688  real(kind=dp) function string_to_real(string)
689 
690  !> String to convert.
691  character(len=*), intent(in) :: string
692 
693  real(kind=dp) :: val
694  integer :: ios
695  character(len=len(string)+300) :: error_msg
696 
697  call assert_msg(733728030, len_trim(string) <= 30, &
698  'error converting "' // trim(string) // '" to real: string too long')
699  read(string, '(f30.0)', iostat=ios) val
700  call assert_msg(727430097, ios == 0, &
701  'error converting "' // trim(string) &
702  // '" to real: IOSTAT = ' // trim(integer_to_string(ios)))
703  string_to_real = val
704 
705  end function string_to_real
706 
707 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
708 
709  !> Convert a string to a logical.
710  logical function string_to_logical(string)
711 
712  !> String to convert.
713  character(len=*), intent(in) :: string
714 
715  logical :: val
716  integer :: ios
717  character(len=len(string)+300) :: error_msg
718 
719  val = .false.
720  if ((trim(string) == 'yes') &
721  .or. (trim(string) == 'y') &
722  .or. (trim(string) == 'true') &
723  .or. (trim(string) == 't') &
724  .or. (trim(string) == '1')) then
725  val = .true.
726  elseif ((trim(string) == 'no') &
727  .or. (trim(string) == 'n') &
728  .or. (trim(string) == 'false') &
729  .or. (trim(string) == 'f') &
730  .or. (trim(string) == '0')) then
731  val = .false.
732  else
733  call die_msg(985010153, 'error converting "' // trim(string) &
734  // '" to logical')
735  end if
736  string_to_logical = val
737 
738  end function string_to_logical
739 
740 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
741 
742  !> Convert an integer to a string format.
743  character(len=PMC_UTIL_CONVERT_STRING_LEN) function integer_to_string(val)
744 
745  !> Value to convert.
746  integer, intent(in) :: val
747 
748  character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
749 
750  ret_val = ""
751  write(ret_val, '(i30)') val
752  integer_to_string = adjustl(ret_val)
753 
754  end function integer_to_string
755 
756 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
757 
758  !> Convert a real to a string format.
759  character(len=PMC_UTIL_CONVERT_STRING_LEN) function real_to_string(val)
760 
761  !> Value to convert.
762  real(kind=dp), intent(in) :: val
763 
764  character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
765 
766  ret_val = ""
767  write(ret_val, '(g30.20)') val
768  real_to_string = adjustl(ret_val)
769 
770  end function real_to_string
771 
772 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
773 
774  !> Convert a logical to a string format.
775  character(len=PMC_UTIL_CONVERT_STRING_LEN) function logical_to_string(val)
776 
777  !> Value to convert.
778  logical, intent(in) :: val
779 
780  character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
781 
782  ret_val = ""
783  if (val) then
784  ret_val = "TRUE"
785  else
786  ret_val = "FALSE"
787  end if
788  logical_to_string = ret_val
789 
790  end function logical_to_string
791 
792 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
793 
794  !> Convert a complex to a string format.
795  character(len=PMC_UTIL_CONVERT_STRING_LEN) function complex_to_string(val)
796 
797  !> Value to convert.
798  complex(kind=dc), intent(in) :: val
799 
800  character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
801 
802  ret_val = ""
803  ret_val = "(" // trim(real_to_string(real(val))) &
804  // ", " // trim(real_to_string(aimag(val))) // ")"
805  complex_to_string = adjustl(ret_val)
806 
807  end function complex_to_string
808 
809 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
810 
811  !> Convert an integer to a string format of maximum length.
812  character(len=PMC_UTIL_CONVERT_STRING_LEN) &
813  function integer_to_string_max_len(val, max_len)
814 
815  !> Value to convert.
816  integer, intent(in) :: val
817  !> Maximum length of resulting string.
818  integer, intent(in) :: max_len
819 
820  character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
821 
822  ret_val = integer_to_string(val)
823  if (len_trim(ret_val) > max_len) then
824  ret_val = real_to_string_max_len(real(val, kind=dp), max_len)
825  end if
826  integer_to_string_max_len = adjustl(ret_val)
827 
828  end function integer_to_string_max_len
829 
830 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
831 
832  !> Convert a real to a string format of maximum length.
833  character(len=PMC_UTIL_CONVERT_STRING_LEN) &
834  function real_to_string_max_len(val, max_len)
835 
836  !> Value to convert.
837  real(kind=dp), intent(in) :: val
838  !> Maximum length of resulting string.
839  integer, intent(in) :: max_len
840 
841  character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val, exp_str, frac_str
842  integer :: exp_val, exp_len, frac_len, use_frac_len, min_frac_len, i
843  real(kind=dp) :: frac_val
844 
845  ret_val = ""
846  if (val == 0d0) then
847  if (max_len >= 3) then
848  ret_val = "0e0"
849  else
850  do i = 1,max_len
851  ret_val(i:i) = "*"
852  end do
853  end if
854  real_to_string_max_len = adjustl(ret_val)
855  return
856  end if
857 
858  exp_val = floor(log10(abs(val)))
859  frac_val = val / 10d0**exp_val
860  exp_str = integer_to_string(exp_val)
861  frac_str = real_to_string(frac_val)
862 
863  exp_len = len_trim(exp_str)
864  frac_len = len_trim(frac_str)
865  use_frac_len = max_len - 1 - exp_len
866  if (use_frac_len > frac_len) then
867  use_frac_len = frac_len
868  end if
869  if (val < 0d0) then
870  min_frac_len = 2
871  else
872  min_frac_len = 1
873  end if
874  if (use_frac_len < min_frac_len) then
875  do i = 1,max_len
876  ret_val(i:i) = "*"
877  end do
878  else
879  ret_val = frac_str(1:use_frac_len) // "e" // trim(exp_str)
880  end if
881  real_to_string_max_len = adjustl(ret_val)
882 
883  end function real_to_string_max_len
884 
885 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
886 
887  !> Convert a time to a string format of maximum length.
888  character(len=PMC_UTIL_CONVERT_STRING_LEN) &
889  function time_to_string_max_len(time, max_len)
890 
891  !> Time to convert (s).
892  real(kind=dp), intent(in) :: time
893  !> Maximum length of resulting string.
894  integer, intent(in) :: max_len
895 
896  integer, dimension(4), parameter :: scale = (/ 1, 60, 60, 24 /)
897  character, dimension(4), parameter :: unit = (/ "s", "m", "h", "d" /)
898 
899  character(len=PMC_UTIL_CONVERT_STRING_LEN) :: ret_val
900  integer :: i
901  logical :: len_ok
902  real(kind=dp) :: scaled_time
903 
904  scaled_time = time
905  len_ok = .false.
906  do i = 1,4
907  scaled_time = scaled_time / real(scale(i), kind=dp)
908  ret_val = trim(integer_to_string(nint(scaled_time))) // unit(i)
909  if (len_trim(ret_val) <= max_len) then
910  len_ok = .true.
911  exit
912  end if
913  end do
914  if (.not. len_ok) then
915  ret_val = ""
916  do i = 1,max_len
917  ret_val(i:i) = "*"
918  end do
919  end if
920  time_to_string_max_len = adjustl(ret_val)
921 
922  end function time_to_string_max_len
923 
924 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
925 
926  !> Convert a real-valued vector into an integer-valued vector.
927  !!
928  !! Use n_samp samples. Returns discrete vector whose relative entry
929  !! sizes are \f$ \ell_1 \f$ closest to the continuous vector.
930  subroutine vec_cts_to_disc(n, vec_cts, n_samp, vec_disc)
931 
932  !> Number of entries in vectors.
933  integer, intent(in) :: n
934  !> Continuous vector.
935  real(kind=dp), intent(in) :: vec_cts(n)
936  !> Number of discrete samples to use.
937  integer, intent(in) :: n_samp
938  !> Discretized vector.
939  integer, intent(out) :: vec_disc(n)
940 
941  integer :: k(1)
942  real(kind=dp) :: vec_tot
943 
944  vec_tot = sum(vec_cts)
945 
946  ! assign a best guess for each bin independently
947  vec_disc = nint(vec_cts / vec_tot * real(n_samp, kind=dp))
948 
949  ! if we have too few particles then add more
950  do while (sum(vec_disc) < n_samp)
951  k = minloc(abs(real(vec_disc + 1, kind=dp) - vec_cts) &
952  - abs(real(vec_disc, kind=dp) - vec_cts))
953  vec_disc(k) = vec_disc(k) + 1
954  end do
955 
956  ! if we have too many particles then remove some
957  do while (sum(vec_disc) > n_samp)
958  k = minloc(abs(real(vec_disc - 1, kind=dp) - vec_cts) &
959  - abs(real(vec_disc, kind=dp) - vec_cts))
960  vec_disc(k) = vec_disc(k) - 1
961  end do
962 
963  call assert_msg(323412496, sum(vec_disc) == n_samp, &
964  'generated incorrect number of samples')
965 
966  end subroutine vec_cts_to_disc
967 
968 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
969 
970  !> Computes the average of an array of integer numbers.
971  subroutine average_integer(int_vec, int_avg)
972 
973  !> Array of integer numbers.
974  integer, intent(in) :: int_vec(:)
975  !> Average of int_vec.
976  integer, intent(out) :: int_avg
977 
978  int_avg = sum(int_vec) / size(int_vec)
979 
980  end subroutine average_integer
981 
982 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
983 
984  !> Computes the average of an array of real numbers.
985  subroutine average_real(real_vec, real_avg)
986 
987  !> Array of real numbers.
988  real(kind=dp), intent(in) :: real_vec(:)
989  !> Average of real_vec.
990  real(kind=dp), intent(out) :: real_avg
991 
992  real_avg = sum(real_vec) / real(size(real_vec), kind=dp)
993 
994  end subroutine average_real
995 
996 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
997 
998  !> Allocate or reallocate the given array to ensure it is of the
999  !> given size, preserving any data and/or initializing to 0.
1000  subroutine ensure_real_array_size(x, n, only_grow)
1001 
1002  !> Array of real numbers.
1003  real(kind=dp), intent(inout), allocatable :: x(:)
1004  !> Desired size of array.
1005  integer, intent(in) :: n
1006  !> Whether to only increase the array size (default .true.).
1007  logical, intent(in), optional :: only_grow
1008 
1009  integer :: new_n
1010  real(kind=dp), allocatable :: tmp_x(:)
1011 
1012  if (allocated(x)) then
1013  new_n = n
1014  if (present(only_grow)) then
1015  new_n = max(new_n, size(x))
1016  end if
1017  if (size(x) /= new_n) then
1018  allocate(tmp_x(new_n))
1019  tmp_x = 0d0
1020  tmp_x(1:min(new_n, size(x))) = x(1:min(new_n, size(x)))
1021  call move_alloc(tmp_x, x)
1022  end if
1023  else
1024  allocate(x(n))
1025  x = 0d0
1026  end if
1027 
1028  end subroutine ensure_real_array_size
1029 
1030 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1031 
1032  !> Allocate or reallocate the given array to ensure it is of the
1033  !> given size, preserving any data and/or initializing to 0.
1034  subroutine ensure_real_array_2d_size(x, n1, n2, only_grow)
1035 
1036  !> Array of real numbers.
1037  real(kind=dp), intent(inout), allocatable :: x(:, :)
1038  !> Desired first size of array.
1039  integer, intent(in) :: n1
1040  !> Desired second size of array.
1041  integer, intent(in) :: n2
1042  !> Whether to only increase the array size (default .true.).
1043  logical, intent(in), optional :: only_grow
1044 
1045  integer :: new_n1, new_n2, n1_min, n2_min
1046  real(kind=dp), allocatable :: tmp_x(:, :)
1047 
1048  if (allocated(x)) then
1049  new_n1 = n1
1050  new_n2 = n2
1051  if (present(only_grow)) then
1052  new_n1 = max(new_n1, size(x, 1))
1053  new_n2 = max(new_n2, size(x, 2))
1054  end if
1055  if ((size(x, 1) /= new_n1) .or. (size(x, 2) /= new_n2)) then
1056  allocate(tmp_x(new_n1, new_n2))
1057  n1_min = min(new_n1, size(x, 1))
1058  n2_min = min(new_n2, size(x, 2))
1059  tmp_x = 0d0
1060  tmp_x(1:n1_min, 1:n2_min) = x(1:n1_min, 1:n2_min)
1061  call move_alloc(tmp_x, x)
1062  end if
1063  else
1064  allocate(x(n1, n2))
1065  x = 0d0
1066  end if
1067 
1068  end subroutine ensure_real_array_2d_size
1069 
1070 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1071 
1072  !> Allocate or reallocate the given array to ensure it is of the
1073  !> given size, preserving any data and/or initializing to 0.
1074  subroutine ensure_integer_array_size(x, n, only_grow)
1075 
1076  !> Array of integer numbers.
1077  integer, intent(inout), allocatable :: x(:)
1078  !> Desired size of array.
1079  integer, intent(in) :: n
1080  !> Whether to only increase the array size (default .true.).
1081  logical, intent(in), optional :: only_grow
1082 
1083  integer :: new_n
1084  integer, allocatable :: tmp_x(:)
1085 
1086  if (allocated(x)) then
1087  new_n = n
1088  if (present(only_grow)) then
1089  new_n = max(new_n, size(x))
1090  end if
1091  if (size(x) /= new_n) then
1092  allocate(tmp_x(new_n))
1093  tmp_x = 0
1094  tmp_x(1:min(new_n, size(x))) = x(1:min(new_n, size(x)))
1095  call move_alloc(tmp_x, x)
1096  end if
1097  else
1098  allocate(x(n))
1099  x = 0
1100  end if
1101 
1102  end subroutine ensure_integer_array_size
1103 
1104 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1105 
1106  !> Allocate or reallocate the given array to ensure it is of the
1107  !> given size, preserving any data and/or initializing to 0.
1108  subroutine ensure_integer_array_2d_size(x, n1, n2, only_grow)
1109 
1110  !> Array of integer numbers.
1111  integer, intent(inout), allocatable :: x(:, :)
1112  !> Desired first size of array.
1113  integer, intent(in) :: n1
1114  !> Desired second size of array.
1115  integer, intent(in) :: n2
1116  !> Whether to only increase the array size (default .true.).
1117  logical, intent(in), optional :: only_grow
1118 
1119  integer :: new_n1, new_n2, n1_min, n2_min
1120  integer, allocatable :: tmp_x(:, :)
1121 
1122  if (allocated(x)) then
1123  new_n1 = n1
1124  new_n2 = n2
1125  if (present(only_grow)) then
1126  new_n1 = max(new_n1, size(x, 1))
1127  new_n2 = max(new_n2, size(x, 2))
1128  end if
1129  if ((size(x, 1) /= new_n1) .or. (size(x, 2) /= new_n2)) then
1130  allocate(tmp_x(new_n1, new_n2))
1131  n1_min = min(new_n1, size(x, 1))
1132  n2_min = min(new_n2, size(x, 2))
1133  tmp_x = 0
1134  tmp_x(1:n1_min, 1:n2_min) = x(1:n1_min, 1:n2_min)
1135  call move_alloc(tmp_x, x)
1136  end if
1137  else
1138  allocate(x(n1, n2))
1139  x = 0
1140  end if
1141 
1142  end subroutine ensure_integer_array_2d_size
1143 
1144 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1145 
1146  !> Allocate or reallocate the given array to ensure it is of the
1147  !> given size, without preserving data.
1148  subroutine ensure_string_array_size(x, n)
1149 
1150  !> Array of strings numbers.
1151  character(len=*), intent(inout), allocatable :: x(:)
1152  !> Desired size of array.
1153  integer, intent(in) :: n
1154 
1155  if (allocated(x)) then
1156  if (size(x) /= n) then
1157  deallocate(x)
1158  allocate(x(n))
1159  end if
1160  else
1161  allocate(x(n))
1162  end if
1163 
1164  end subroutine ensure_string_array_size
1165 
1166 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1167 
1168  !> Strip the extension to find the basename.
1169  !!
1170  !! The filename is assumed to be of the form basename.extension,
1171  !! where extension contains no periods. If there is no period in
1172  !! filename then basename is the whole filename.
1173  subroutine get_basename(filename, basename)
1174 
1175  !> Full filename.
1176  character(len=*), intent(in) :: filename
1177  !> Basename.
1178  character(len=*), intent(out) :: basename
1179 
1180  integer :: i
1181  logical :: found_period
1182 
1183  basename = filename
1184  i = len_trim(basename)
1185  found_period = .false.
1186  do while ((i > 0) .and. (.not. found_period))
1187  ! Fortran .and. does not short-circuit, so we can't do the
1188  ! obvious do while ((i > 0) .and. (basename(i:i) /= ".")),
1189  ! instead we have to use this hack with the found_period
1190  ! logical variable.
1191  if (basename(i:i) == ".") then
1192  found_period = .true.
1193  else
1194  i = i - 1
1195  end if
1196  end do
1197  if (i > 0) then
1198  basename(i:) = ""
1199  end if
1200 
1201  end subroutine get_basename
1202 
1203 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1204 
1205  !> Current date and time in ISO 8601 format.
1206  subroutine iso8601_date_and_time(date_time)
1207 
1208  !> Result string.
1209  character(len=*), intent(out) :: date_time
1210 
1211  character(len=10) :: date, time, zone
1212 
1213  call assert_msg(893219839, len(date_time) >= 29, &
1214  "date_time string must have length at least 29")
1215  call date_and_time(date, time, zone)
1216  date_time = ""
1217  write(date_time, '(14a)') date(1:4), "-", date(5:6), "-", &
1218  date(7:8), "T", time(1:2), ":", time(3:4), ":", &
1219  time(5:10), zone(1:3), ":", zone(4:5)
1220 
1221  end subroutine iso8601_date_and_time
1222 
1223 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1224 
1225  !> Convert degrees to radians.
1226  real(kind=dp) function deg2rad(deg)
1227 
1228  !> Input degrees.
1229  real(kind=dp), intent(in) :: deg
1230 
1231  deg2rad = deg / 180d0 * const%pi
1232 
1233  end function deg2rad
1234 
1235 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1236 
1237  !> Convert radians to degrees.
1238  real(kind=dp) function rad2deg(rad)
1239 
1240  !> Input radians.
1241  real(kind=dp), intent(in) :: rad
1242 
1243  rad2deg = rad / const%pi * 180d0
1244 
1245  end function rad2deg
1246 
1247 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1248 
1249  !> Sort the given data array and return the permutation defining the sort.
1250  !!
1251  !! If \c data is the original data array, \c new_data is the sorted
1252  !! value of data, and \c perm is the returned permutation, then
1253  !! <tt>new_data(i) = data(perm(i))</tt>.
1254  subroutine integer_sort(data, perm)
1255 
1256  !> Data array to sort, sorted on exit.
1257  integer, intent(inout) :: data(:)
1258  !> Permutation defining the sort: <tt>new_data(i) = data(perm(i))</tt>.
1259  integer, intent(out) :: perm(size(data))
1260 
1261 #ifdef PMC_USE_C_SORT
1262  integer(kind=c_int) :: n_c
1263  integer(kind=c_int), target :: data_c(size(data))
1264  integer(kind=c_int), target :: perm_c(size(data))
1265  type(c_ptr) :: data_ptr, perm_ptr
1266 
1267 #ifndef DOXYGEN_SKIP_DOC
1268  interface
1269  subroutine integer_sort_c(n_c, data_ptr, perm_ptr) bind(c)
1270  use iso_c_binding
1271  integer(kind=c_int), value :: n_c
1272  type(c_ptr), value :: data_ptr, perm_ptr
1273  end subroutine integer_sort_c
1274  end interface
1275 #endif
1276 
1277  data_c = int(data, kind=c_int)
1278  perm_c = 0_c_int
1279  n_c = int(size(data), kind=c_int)
1280  data_ptr = c_loc(data_c)
1281  perm_ptr = c_loc(perm_c)
1282  call integer_sort_c(n_c, data_ptr, perm_ptr)
1283  data = int(data_c)
1284  perm = int(perm_c)
1285 #endif
1286 
1287  end subroutine integer_sort
1288 
1289 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1290 
1291  !> Load a real array from a text file.
1292  subroutine loadtxt(filename, data)
1293 
1294  !> Filename to read from.
1295  character(len=*), intent(in) :: filename
1296  !> Array to store data in.
1297  real(kind=dp), intent(inout), allocatable :: data(:,:)
1298 
1299  integer :: unit, row, col
1300  logical :: eol, eof
1301  character(len=1000) :: word
1302 
1303  deallocate(data)
1304  allocate(data(1,0))
1305  call open_file_read(filename, unit)
1306 
1307  eof = .false.
1308  row = 1
1309  col = 1
1310  do while (.not. eof)
1311  call read_word_raw(unit, word, eol, eof)
1312  if (len_trim(word) > 0) then
1313  if (row == 1) then
1314  if (col > size(data, 2)) then
1315  call reallocate_real_array2d(data, 1, 2 * col)
1316  end if
1317  else
1318  if (col > size(data, 2)) then
1319  call assert_msg(516120334, col <= size(data, 2), &
1320  trim(filename) // ": line " &
1321  // trim(integer_to_string(row)) // " too long")
1322  end if
1323  end if
1324  if (row > size(data, 1)) then
1325  call reallocate_real_array2d(data, 2 * row, size(data, 2))
1326  end if
1327  data(row, col) = string_to_real(word)
1328  col = col + 1
1329  end if
1330  if (eol .or. eof) then
1331  if (row == 1) then
1332  call reallocate_real_array2d(data, 1, col - 1)
1333  else
1334  call assert_msg(266924956, &
1335  (col == 1) .or. (col == size(data, 2) + 1), &
1336  trim(filename) // ": line " &
1337  // trim(integer_to_string(row)) // " too short")
1338  end if
1339  end if
1340  if (eol) then
1341  row = row + 1
1342  col = 1
1343  end if
1344  end do
1345 
1346  if (col == 1) then
1347  call reallocate_real_array2d(data, row - 1, size(data, 2))
1348  end if
1349  call close_file(unit)
1350 
1351  end subroutine loadtxt
1352 
1353 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1354 
1355  !> Reallocate a 2D real array to the given size, preserving the contents.
1356  subroutine reallocate_real_array2d(data, rows, cols)
1357 
1358  !> Array to reallocate.
1359  real(kind=dp), allocatable, intent(inout) :: data(:,:)
1360  !> New leading dimension.
1361  integer, intent(in) :: rows
1362  !> New trailing dimension.
1363  integer, intent(in) :: cols
1364 
1365  real(kind=dp) :: tmp_data(rows, cols)
1366  integer :: data_rows, data_cols
1367 
1368  data_rows = min(rows, size(data, 1))
1369  data_cols = min(cols, size(data, 2))
1370  tmp_data(1:data_rows, 1:data_cols) = data(1:data_rows, 1:data_cols)
1371  deallocate(data)
1372  allocate(data(rows, cols))
1373  data(1:data_rows, 1:data_cols) = tmp_data(1:data_rows, 1:data_cols)
1374 
1375  end subroutine reallocate_real_array2d
1376 
1377 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1378 
1379  !> Read a single character from a file, signaling if we have hit
1380  !> end-of-line (EOL) or end-of-file (EOF). If EOL or EOF are true
1381  !> then the character value should be ignored.
1382  !!
1383  !! Testing with gfortran 4.5.3 and different length files shows:
1384  !!
1385  !! Empty file (total file length 0):
1386  !! * eol = .false., eof = .true.
1387  !! * subsequent calls error
1388  !!
1389  !! File containing a single 'A' character (total file length 1):
1390  !! * char = 'A', eol = .false., eof = .false.
1391  !! * eol = .false., eof = .true.
1392  !! * subsequent calls error
1393  !!
1394  !! File containing a single newline '\n' (total file length 1):
1395  !! * eol = .true., eof = .false.
1396  !! * eol = .false., eof = .true.
1397  !! * subsequent calls error
1398  !!
1399  !! File containing a character and newline 'A\n' (total file length 2):
1400  !! * char = 'A', eol = .false., eof = .false.
1401  !! * eol = .true., eof = .false.
1402  !! * eol = .false., eof = .true.
1403  !! * subsequent calls error
1404  subroutine read_char_raw(unit, char, eol, eof)
1405 
1406  !> Unit number to read from.
1407  integer, intent(in) :: unit
1408  !> Character read.
1409  character, intent(out) :: char
1410  !> True if at EOL (end of line).
1411  logical, intent(out) :: eol
1412  !> True if at EOF (end of file).
1413  logical, intent(out) :: eof
1414 
1415  integer :: ios, n_read
1416  character(len=1) :: read_char
1417 
1418  eol = .false.
1419  eof = .false.
1420  char = " " ! shut up uninitialized variable warnings
1421  read_char = "" ! needed for pgf95 for reading blank lines
1422  read(unit=unit, fmt='(a)', advance='no', end=100, eor=110, &
1423  iostat=ios) read_char
1424  if (ios /= 0) then
1425  write(0,*) 'ERROR: reading file: IOSTAT = ', ios
1426  stop 2
1427  end if
1428  ! only reach here if we didn't hit end-of-record (end-of-line) in
1429  ! the above read
1430  char = read_char
1431  goto 120
1432 
1433 100 eof = .true. ! goto here if end-of-file was encountered immediately
1434  goto 120
1435 
1436 110 eol = .true. ! goto here if end-of-record, meaning end-of-line
1437 
1438 120 return
1439 
1440  end subroutine read_char_raw
1441 
1442 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1443 
1444  !> Read a white-space delimited word from a file, signaling if we
1445  !> have EOL or EOF. If EOL or EOF are true then the word will still
1446  !> be meaningful data. If there was no data to be read then
1447  !> len(word) will be 0.
1448  subroutine read_word_raw(unit, word, eol, eof)
1449 
1450  !> Unit number to read from.
1451  integer, intent(in) :: unit
1452  !> Word read.
1453  character(len=*), intent(out) :: word
1454  !> True if at EOL (end of line).
1455  logical, intent(out) :: eol
1456  !> True if at EOF (end of file).
1457  logical, intent(out) :: eof
1458 
1459  integer :: i
1460  character :: char
1461 
1462  word = ""
1463 
1464  ! skip over spaces
1465  call read_char_raw(unit, char, eol, eof)
1466  do while (((ichar(char) == 9) .or. (ichar(char) == 32)) &
1467  .and. (.not. eol) .and. (.not. eof))
1468  call read_char_raw(unit, char, eol, eof)
1469  end do
1470  if (eol .or. eof) return
1471 
1472  ! char is now the first word character
1473  i = 1
1474  word(i:i) = char
1475  call read_char_raw(unit, char, eol, eof)
1476  do while ((ichar(char) /= 9) .and. (ichar(char) /= 32) &
1477  .and. (.not. eol) .and. (.not. eof) .and. (i < len(word)))
1478  i = i + 1
1479  word(i:i) = char
1480  if (i < len(word)) then
1481  call read_char_raw(unit, char, eol, eof)
1482  end if
1483  end do
1484 
1485  end subroutine read_word_raw
1486 
1487 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1488 
1489  !> Checks whether a string starts with a given other string.
1490  !!
1491  !! <tt>starts_with(A, B)</tt> returns \c true if string \c A starts
1492  !! with string \c B.
1493  logical function starts_with(string, start_string)
1494 
1495  !> String to test.
1496  character(len=*), intent(in) :: string
1497  !> Starting string.
1498  character(len=*), intent(in) :: start_string
1499 
1500  if (len(string) < len(start_string)) then
1501  starts_with = .false.
1502  return
1503  end if
1504  if (string(1:len(start_string)) == start_string) then
1505  starts_with = .true.
1506  else
1507  starts_with = .false.
1508  end if
1509 
1510  end function starts_with
1511 
1512 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1513 
1514  !> Copy a 1D array of reals, reallocating if necessary.
1515  subroutine copy_real_1d(source, dest)
1516 
1517  !> Source array.
1518  real(kind=dp), intent(in) :: source(:)
1519  !> Destination array.
1520  real(kind=dp), intent(inout), pointer :: dest(:)
1521 
1522  if (size(dest) /= size(source)) then
1523  deallocate(dest)
1524  allocate(dest(size(source)))
1525  end if
1526  dest = source
1527 
1528  end subroutine copy_real_1d
1529 
1530 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1531 
1532  !> Copy a 2D array of reals, reallocating if necessary.
1533  subroutine copy_real_2d(source, dest)
1534 
1535  !> Source array.
1536  real(kind=dp), intent(in) :: source(:, :)
1537  !> Destination array.
1538  real(kind=dp), intent(inout), pointer :: dest(:, :)
1539 
1540  if ((size(dest, 1) /= size(source, 1)) &
1541  .or. (size(dest, 2) /= size(source, 2))) then
1542  deallocate(dest)
1543  allocate(dest(size(source, 1), size(source, 2)))
1544  end if
1545  dest = source
1546 
1547  end subroutine copy_real_2d
1548 
1549 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1550 
1551  !> Compute the harmonic mean of two numbers.
1552  elemental real(kind=dp) function harmonic_mean(x1, x2)
1553 
1554  !> First number to average.
1555  real(kind=dp), intent(in) :: x1
1556  !> Second number to average.
1557  real(kind=dp), intent(in) :: x2
1558 
1559  harmonic_mean = 1d0 / (1d0 / x1 + 1d0 / x2)
1560 
1561  end function harmonic_mean
1562 
1563 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1564 
1565  !> Compute \f$ - p \ln p\f$ for computing entropy.
1566  elemental real(kind=dp) function nplogp(p)
1567 
1568  !> Probability \f$p\f$.
1569  real(kind=dp), intent(in) :: p
1570 
1571  if (p <= 0d0) then
1572  nplogp = 0d0
1573  else
1574  nplogp = - p * log(p)
1575  end if
1576 
1577  end function nplogp
1578 
1579 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1580 
1581  !> Compute the entropy of a probability mass function (non
1582  !> necessarily normalized).
1583  real(kind=dp) function entropy(p)
1584 
1585  !> Probability mass function, will be normalized before use.
1586  real(kind=dp), intent(in) :: p(:)
1587 
1588  entropy = sum(nplogp(p / sum(p)))
1589 
1590  end function entropy
1591 
1592 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1593 
1594 end module pmc_util
real(kind=dp) elemental function vol2diam(v)
Convert volume (m^3) to diameter (m).
Definition: util.F90:250
character(len=pmc_util_convert_string_len) function integer_to_string_max_len(val, max_len)
Convert an integer to a string format of maximum length.
Definition: util.F90:812
subroutine read_char_raw(unit, char, eol, eof)
Read a single character from a file, signaling if we have hit EOL or EOF. If EOL or EOF are true then...
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:133
real(kind=dp) function rad2deg(rad)
Convert radians to degrees.
Definition: util.F90:1238
subroutine copy_real_2d(source, dest)
Copy a 2D array of reals, reallocating if necessary.
Definition: util.F90:1533
subroutine ensure_string_array_size(x, n)
Allocate or reallocate the given array to ensure it is of the given size, without preserving data...
Definition: util.F90:1148
integer function linspace_find(min_x, max_x, n, x)
Find the position of a real number in a 1D linear array.
Definition: util.F90:517
real(kind=dp) elemental function diam2vol(d)
Convert diameter (m) to volume (m^3).
Definition: util.F90:298
subroutine iso8601_date_and_time(date_time)
Current date and time in ISO 8601 format.
Definition: util.F90:1206
real(kind=dp) elemental function rad2vol(r)
Convert radius (m) to volume (m^3).
Definition: util.F90:274
real(kind=dp) function entropy(p)
Compute the entropy of a probability mass function (non necessarily normalized).
Definition: util.F90:1583
int integer_sort_c(int n, int *data, int *perm)
Sort the given data array and return the permutation.
Definition: sort.c:40
character(len=pmc_util_convert_string_len) function logical_to_string(val)
Convert a logical to a string format.
Definition: util.F90:775
logical function starts_with(string, start_string)
Checks whether a string starts with a given other string.
Definition: util.F90:1493
Physical constants.
Definition: constants.F90:9
subroutine close_file(unit)
Close a file and de-assign the unit.
Definition: util.F90:225
integer function find_1d(n, x_vals, x)
Find the position of a real number in an arbitrary 1D array.
Definition: util.F90:577
real(kind=dp) elemental function rad2diam(r)
Convert radius (m) to diameter (m).
Definition: util.F90:262
subroutine ensure_integer_array_2d_size(x, n1, n2, only_grow)
Allocate or reallocate the given array to ensure it is of the given size, preserving any data and/or ...
Definition: util.F90:1108
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:76
subroutine logspace(min_x, max_x, x)
Makes a logarithmically spaced array of length n from min to max.
Definition: util.F90:481
logical function almost_equal(d1, d2)
Tests whether two real numbers are almost equal using only a relative tolerance.
Definition: util.F90:311
real(kind=dp) function interp_linear_disc(x_1, x_n, n, i)
Linear interpolation over discrete indices.
Definition: util.F90:644
real(kind=dp) elemental function diam2rad(d)
Convert diameter (m) to radius (m).
Definition: util.F90:286
integer function string_to_integer(string)
Convert a string to an integer.
Definition: util.F90:665
subroutine check_time_multiple(first_name, first_time, second_name, second_time)
Check that the first time interval is close to an integer multiple of the second, and warn if it is n...
Definition: util.F90:368
subroutine average_integer(int_vec, int_avg)
Computes the average of an array of integer numbers.
Definition: util.F90:971
subroutine warn_assert_msg(code, condition_ok, warning_msg)
Prints a warning message if condition_ok is false.
Definition: util.F90:58
subroutine copy_real_1d(source, dest)
Copy a 1D array of reals, reallocating if necessary.
Definition: util.F90:1515
Common utility subroutines.
Definition: util.F90:9
subroutine reallocate_real_array2d(data, rows, cols)
Reallocate a 2D real array to the given size, preserving the contents.
Definition: util.F90:1356
logical function string_to_logical(string)
Convert a string to a logical.
Definition: util.F90:710
subroutine average_real(real_vec, real_avg)
Computes the average of an array of real numbers.
Definition: util.F90:985
real(kind=dp) function deg2rad(deg)
Convert degrees to radians.
Definition: util.F90:1226
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:743
character(len=pmc_util_convert_string_len) function complex_to_string(val)
Convert a complex to a string format.
Definition: util.F90:795
logical function almost_equal_abs(d1, d2, abs_tol)
Tests whether two real numbers are almost equal using an absolute and relative tolerance.
Definition: util.F90:338
subroutine die(code)
Error immediately.
Definition: util.F90:121
integer function get_unit()
Returns an available unit number. This should be freed by free_unit().
Definition: util.F90:147
subroutine ensure_real_array_size(x, n, only_grow)
Allocate or reallocate the given array to ensure it is of the given size, preserving any data and/or ...
Definition: util.F90:1000
subroutine warn_msg(code, warning_msg, already_warned)
Prints a warning message.
Definition: util.F90:36
subroutine loadtxt(filename, data)
Load a real array from a text file.
Definition: util.F90:1292
subroutine check_event(time, timestep, interval, last_time, do_event)
Computes whether an event is scheduled to take place.
Definition: util.F90:399
subroutine ensure_real_array_2d_size(x, n1, n2, only_grow)
Allocate or reallocate the given array to ensure it is of the given size, preserving any data and/or ...
Definition: util.F90:1034
elemental real(kind=dp) function harmonic_mean(x1, x2)
Compute the harmonic mean of two numbers.
Definition: util.F90:1552
subroutine integer_sort(data, perm)
Sort the given data array and return the permutation defining the sort.
Definition: util.F90:1254
real(kind=dp) function interp_1d(x_vals, y_vals, x)
1D linear interpolation.
Definition: util.F90:611
real(kind=dp) elemental function vol2rad(v)
Convert volume (m^3) to radius (m).
Definition: util.F90:238
character(len=pmc_util_convert_string_len) function real_to_string(val)
Convert a real to a string format.
Definition: util.F90:759
integer function logspace_find(min_x, max_x, n, x)
Find the position of a real number in a 1D logarithmic array.
Definition: util.F90:551
subroutine linspace(min_x, max_x, x)
Makes a linearly spaced array from min to max.
Definition: util.F90:453
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:102
subroutine free_unit(unit)
Frees a unit number returned by get_unit().
Definition: util.F90:171
real(kind=dp) function string_to_real(string)
Convert a string to a real.
subroutine read_word_raw(unit, word, eol, eof)
Read a white-space delimited word from a file, signaling if we have EOL or EOF. If EOL or EOF are tru...
character(len=pmc_util_convert_string_len) function time_to_string_max_len(time, max_len)
Convert a time to a string format of maximum length.
Definition: util.F90:888
subroutine vec_cts_to_disc(n, vec_cts, n_samp, vec_disc)
Convert a real-valued vector into an integer-valued vector.
Definition: util.F90:930
character(len=pmc_util_convert_string_len) function real_to_string_max_len(val, max_len)
Convert a real to a string format of maximum length.
Definition: util.F90:833
subroutine open_file_read(filename, unit)
Open a file for reading with an automatically assigned unit and test that it succeeds. The file should be closed with close_file().
Definition: util.F90:184
elemental real(kind=dp) function nplogp(p)
Compute for computing entropy.
Definition: util.F90:1566
subroutine ensure_integer_array_size(x, n, only_grow)
Allocate or reallocate the given array to ensure it is of the given size, preserving any data and/or ...
Definition: util.F90:1074
subroutine open_file_write(filename, unit)
Open a file for writing with an automatically assigned unit and test that it succeeds. The file should be closed with close_file().
Definition: util.F90:205
subroutine get_basename(filename, basename)
Strip the extension to find the basename.
Definition: util.F90:1173