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