gray/src/utils.f90
Michele Guerini Rocco d5c81268de
src/utils.f90: clean up
- Replace the `get_free_unit` subroutine with the built-in
  `newutin` option of the `open` statement.

- Replace `locatex` with just `locate` + an index offset.

- Replace `inside` with `contour%contains`.

- Merge `vmaxmin` and `vmaxmini` into a single subroutine
  with optional arguments.

- Remove unused `range2rect`, `bubble`.
2024-11-03 09:19:21 +01:00

194 lines
4.7 KiB
Fortran

module utils
use const_and_precisions, only : wp_
implicit none
contains
pure subroutine locate(array, value, index)
! Given an `array`, and a `value` returns the `index`
! such that `array(index) < value < array(index+1)`
!
! Notes:
! 1. `array` must be monotonic, either increasing or decreasing.
! 2. an index equal to 0 or `size(array)` indicate the value was not found
! Source: Numerical Recipes
! subroutine arguments
real(wp_), intent(in) :: array(:), value
integer, intent(out) :: index
! local variables
integer :: jl, ju, jm
logical :: incr
jl = 0
ju = size(array) + 1
incr = array(size(array)) > array(1)
do while ((ju - jl) > 1)
jm = (ju + jl)/2
if (incr .eqv. (value > array(jm))) then
jl = jm
else
ju = jm
endif
end do
index = jl
end subroutine locate
pure subroutine locate_unord(array, value, locs, n, nlocs)
! Given an `array` of size `n` and a `value`, finds at most
! `n` locations `locs` such that `value` is between
! `array(locs(i))` and `array(locs(i+i))`, in whichever order.
! subroutine arguments
real(wp_), intent(in) :: array(:)
real(wp_), intent(in) :: value
integer, intent(inout) :: locs(n)
integer, intent(in) :: n
integer, intent(out) :: nlocs
! local variables
integer :: i
logical :: larger_than_last
nlocs = 0
if (size(array) < 2) return
larger_than_last = value > array(1)
do i = 2, size(array)
! Note: the condition is equivalent to
! (array(i-1) < value < array(i))
! .or. (array(i) < value < array(i-1))
if (array(i) < value .neqv. larger_than_last) then
larger_than_last = value > array(i)
nlocs = nlocs + 1
if (nlocs <= n) locs(nlocs) = i - 1
end if
end do
end subroutine locate_unord
pure function linear_interp(x1, y1, x2, y2, x) result(y)
! Computes the linear interpolation between
! the points `(x1, y1)` and `(x2, y2)` at `x`.
!
! Note: x1 != x2 is assumed
! function arguments
real(wp_), intent(in) :: x1, y1, x2, y2, x
real(wp_) :: y
! local variables
real(wp_) :: a
a = (x2 - x)/(x2 - x1)
y = a*y1 + (1 - a)*y2
end function linear_interp
pure subroutine vmaxmin(x, xmin, xmax, imin, imax)
! Computes the maximum and minimum of the array `x`
! and optionally their indices
! subroutine arguments
real(wp_), intent(in) :: x(:)
real(wp_), intent(out) :: xmin, xmax
integer, intent(out), optional :: imin, imax
! local variables
integer :: i
if (size(x) < 1) then
if (present(imin)) imin = 0
if (present(imax)) imax = 0
return
end if
if (present(imin)) imin = 1
if (present(imax)) imax = 1
xmin = x(1)
xmax = x(1)
do i = 2, size(x)
if(x(i) < xmin) then
xmin = x(i)
if (present(imin)) imin = i
else if(x(i) > xmax) then
xmax = x(i)
if (present(imax)) imax = i
end if
end do
end subroutine vmaxmin
pure subroutine sort_pair(p)
! Sorts the pair `p` in ascending order
! subroutine inputs
real(wp_), intent(inout) :: p(2)
! local variables
real(wp_) :: temp
if (p(1) > p(2)) then
temp = p(1)
p(1) = p(2)
p(2) = temp
end if
end subroutine sort_pair
function dirname(filepath) result(directory)
! Get the parent `directory` of `filepath`
! function arguments
character(*), intent(in) :: filepath
character(:), allocatable :: directory
! local variables
character(255) :: cwd
integer :: last_sep
last_sep = scan(filepath, '/', back=.true.)
directory = filepath(1:last_sep)
! append the cwd to relative paths
if (isrelative(filepath)) then
call getcwd(cwd)
directory = trim(cwd) // '/' // directory
end if
end function dirname
function isrelative(filepath)
! Check if `filepath` is a relative or an absolute path
! function arguments
character(*), intent(in) :: filepath
logical :: isrelative
isrelative = (filepath(1:1) /= '/')
end function isrelative
pure function digits(n)
! Returns the number of digits of an integer
! function arguments
integer, intent(in) :: n
integer :: digits
! How it works:
! 1. `bit_size` returns the number of bits + 1 (sign)
! 2. subtracting the number of leading zeros gives the significant bits
! 3. multiplying by log₁₀(2) gives the approximate number of digits
! 4. `ceiling` rounds to the next closest integer
digits = ceiling((bit_size(abs(n)) - leadz(abs(n))) * log10(2.0_wp_))
end function digits
end module utils