diff --git a/src/stdlib_experimental_io.f90 b/src/stdlib_experimental_io.f90 index 9506bf234..0e0ca6bab 100644 --- a/src/stdlib_experimental_io.f90 +++ b/src/stdlib_experimental_io.f90 @@ -1,5 +1,5 @@ module stdlib_experimental_io -use iso_fortran_env, only: sp=>real32, dp=>real64 +use iso_fortran_env, only: sp=>real32, dp=>real64, qp=>real128 implicit none private public :: loadtxt, savetxt @@ -7,22 +7,58 @@ module stdlib_experimental_io interface loadtxt module procedure sloadtxt module procedure dloadtxt + module procedure qloadtxt end interface interface savetxt module procedure ssavetxt module procedure dsavetxt + module procedure qsavetxt end interface contains subroutine sloadtxt(filename, d) +! Loads a 2D array from a text file. +! +! Arguments +! --------- +! +! Filename to load the array from character(len=*), intent(in) :: filename +! The array 'd' will be automatically allocated with the correct dimensions real(sp), allocatable, intent(out) :: d(:,:) -real(dp), allocatable :: tmp(:,:) -call dloadtxt(filename, tmp) -allocate(d(size(tmp,1),size(tmp,2))) -d = real(tmp,sp) +! +! Example +! ------- +! +! real(sp), allocatable :: data(:, :) +! call loadtxt("log.txt", data) ! 'data' will be automatically allocated +! +! Where 'log.txt' contains for example:: +! +! 1 2 3 +! 2 4 6 +! 8 9 10 +! 11 12 13 +! ... +! +integer :: s +integer :: nrow,ncol,i + +open(newunit=s, file=filename, status="old") + +! determine number of columns +ncol = number_of_columns(s) + +! determine number or rows +nrow = number_of_rows_numeric(s) + +allocate(d(nrow, ncol)) +do i = 1, nrow + read(s, *) d(i, :) +end do +close(s) end subroutine subroutine dloadtxt(filename, d) @@ -50,34 +86,59 @@ subroutine dloadtxt(filename, d) ! 11 12 13 ! ... ! -character :: c -integer :: s, ncol, nrow, ios, i -logical :: lastwhite -real(dp) :: r +integer :: s +integer :: nrow,ncol,i open(newunit=s, file=filename, status="old") ! determine number of columns -ncol = 0 -lastwhite = .true. -do - read(s, '(a)', advance='no', iostat=ios) c - if (ios /= 0) exit - if (lastwhite .and. .not. whitechar(c)) ncol = ncol + 1 - lastwhite = whitechar(c) -end do - -rewind(s) +ncol = number_of_columns(s) ! determine number or rows -nrow = 0 -do - read(s, *, iostat=ios) r - if (ios /= 0) exit - nrow = nrow + 1 +nrow = number_of_rows_numeric(s) + +allocate(d(nrow, ncol)) +do i = 1, nrow + read(s, *) d(i, :) end do +close(s) +end subroutine + +subroutine qloadtxt(filename, d) +! Loads a 2D array from a text file. +! +! Arguments +! --------- +! +! Filename to load the array from +character(len=*), intent(in) :: filename +! The array 'd' will be automatically allocated with the correct dimensions +real(qp), allocatable, intent(out) :: d(:,:) +! +! Example +! ------- +! +! real(qp), allocatable :: data(:, :) +! call loadtxt("log.txt", data) ! 'data' will be automatically allocated +! +! Where 'log.txt' contains for example:: +! +! 1 2 3 +! 2 4 6 +! 8 9 10 +! 11 12 13 +! ... +! +integer :: s +integer :: nrow,ncol,i + +open(newunit=s, file=filename, status="old") -rewind(s) +! determine number of columns +ncol = number_of_columns(s) + +! determine number or rows +nrow = number_of_rows_numeric(s) allocate(d(nrow, ncol)) do i = 1, nrow @@ -86,10 +147,28 @@ subroutine dloadtxt(filename, d) close(s) end subroutine + subroutine ssavetxt(filename, d) -character(len=*), intent(in) :: filename -real(sp), intent(in) :: d(:,:) -call dsavetxt(filename, real(d,dp)) +! Saves a 2D array into a textfile. +! +! Arguments +! --------- +! +character(len=*), intent(in) :: filename ! File to save the array to +real(sp), intent(in) :: d(:,:) ! The 2D array to save +! +! Example +! ------- +! +! real(sp) :: data(3, 2) +! call savetxt("log.txt", data) + +integer :: s, i +open(newunit=s, file=filename, status="replace") +do i = 1, size(d, 1) + write(s, *) d(i, :) +end do +close(s) end subroutine subroutine dsavetxt(filename, d) @@ -115,6 +194,69 @@ subroutine dsavetxt(filename, d) close(s) end subroutine +subroutine qsavetxt(filename, d) +! Saves a 2D array into a textfile. +! +! Arguments +! --------- +! +character(len=*), intent(in) :: filename ! File to save the array to +real(qp), intent(in) :: d(:,:) ! The 2D array to save +! +! Example +! ------- +! +! real(dp) :: data(3, 2) +! call savetxt("log.txt", data) + +integer :: s, i +open(newunit=s, file=filename, status="replace") +do i = 1, size(d, 1) + write(s, *) d(i, :) +end do +close(s) +end subroutine + + +integer function number_of_columns(s) + ! determine number of columns + integer,intent(in)::s + + integer :: ios + character :: c + logical :: lastwhite + + rewind(s) + number_of_columns = 0 + lastwhite = .true. + do + read(s, '(a)', advance='no', iostat=ios) c + if (ios /= 0) exit + if (lastwhite .and. .not. whitechar(c)) number_of_columns = number_of_columns + 1 + lastwhite = whitechar(c) + end do + rewind(s) + +end function + +integer function number_of_rows_numeric(s) + ! determine number or rows + integer,intent(in)::s + integer :: ios + + real::r + + rewind(s) + number_of_rows_numeric = 0 + do + read(s, *, iostat=ios) r + if (ios /= 0) exit + number_of_rows_numeric = number_of_rows_numeric + 1 + end do + + rewind(s) + +end function logical function whitechar(char) ! white character ! returns .true. if char is space (32) or tab (9), .false. otherwise diff --git a/src/tests/loadtxt/test_loadtxt.f90 b/src/tests/loadtxt/test_loadtxt.f90 index 490fcea38..1d07c5e13 100644 --- a/src/tests/loadtxt/test_loadtxt.f90 +++ b/src/tests/loadtxt/test_loadtxt.f90 @@ -1,9 +1,15 @@ program test_loadtxt -use iso_fortran_env, only: dp=>real64 +use iso_fortran_env, only: sp=>real32, dp=>real64 ,qp=>real128 use stdlib_experimental_io, only: loadtxt implicit none +real(sp), allocatable :: s(:, :) real(dp), allocatable :: d(:, :) +real(qp), allocatable :: q(:, :) + +call loadtxt("array1.dat", s) +call print_array(s) + call loadtxt("array1.dat", d) call print_array(d) @@ -16,15 +22,34 @@ program test_loadtxt call loadtxt("array4.dat", d) call print_array(d) +call loadtxt("array4.dat", q) +call print_array(q) + contains subroutine print_array(a) -real(dp) :: a(:, :) +class(*),intent(in) :: a(:, :) integer :: i print *, "Array, shape=(", size(a, 1), ",", size(a, 2), ")" -do i = 1, size(a, 1) + + select type(a) + type is(real(sp)) + do i = 1, size(a, 1) print *, a(i, :) -end do + end do + type is(real(dp)) + do i = 1, size(a, 1) + print *, a(i, :) + end do + type is(real(qp)) + do i = 1, size(a, 1) + print *, a(i, :) + end do + class default + write(*,'(a)')'The proposed type is not supported' + error stop + end select + end subroutine end program diff --git a/src/tests/loadtxt/test_savetxt.f90 b/src/tests/loadtxt/test_savetxt.f90 index ca6344b83..b5d991f6f 100644 --- a/src/tests/loadtxt/test_savetxt.f90 +++ b/src/tests/loadtxt/test_savetxt.f90 @@ -1,11 +1,12 @@ program test_loadtxt -use iso_fortran_env, only: sp=>real32, dp=>real64 +use iso_fortran_env, only: sp=>real32, dp=>real64 ,qp=>real128 use stdlib_experimental_io, only: loadtxt, savetxt use stdlib_experimental_error, only: assert implicit none call test_sp() call test_dp() +call test_qp() contains @@ -42,4 +43,20 @@ subroutine test_dp() call assert(all(abs(e-d2) < epsilon(1._dp))) end subroutine + subroutine test_qp() + real(qp) :: d(3, 2), e(2, 3) + real(qp), allocatable :: d2(:, :) + d = reshape([1, 2, 3, 4, 5, 6], [3, 2]) + call savetxt("tmp.dat", d) + call loadtxt("tmp.dat", d2) + call assert(all(shape(d2) == [3, 2])) + call assert(all(abs(d-d2) < epsilon(1._qp))) + + e = reshape([1, 2, 3, 4, 5, 6], [2, 3]) + call savetxt("tmp.dat", e) + call loadtxt("tmp.dat", d2) + call assert(all(shape(d2) == [2, 3])) + call assert(all(abs(e-d2) < epsilon(1._qp))) + end subroutine + end program