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