1
0
mirror of https://github.com/TREX-CoE/trexio.git synced 2024-12-22 20:35:44 +01:00

Fix fortran tests

This commit is contained in:
Anthony Scemama 2021-06-04 18:14:45 +02:00
parent 31bf17c66a
commit ea9ffaea43
4 changed files with 135 additions and 32 deletions

View File

@ -87,7 +87,7 @@ tests_test_f_SOURCES = $(test_trexio_f) tests/test_f.f90
tests_test_f_LDADD = src/libtrexio.la
tests_test_f_LDFLAGS = -no-install
$(test_trexio_f):
$(test_trexio_f): $(trexio_f)
cp $(trexio_f) $(test_trexio_f)
clean-local:

View File

@ -25,11 +25,12 @@ module trexio
implicit none
integer, parameter :: trexio_exit_code = 4
integer, parameter :: trexio_backend = 4
integer, parameter :: TREXIO_HDF5 = 0
integer, parameter :: TREXIO_TEXT = 1
! integer, parameter :: TREXIO_JSON = 2
integer, parameter :: TREXIO_INVALID_BACK_END = 2
integer(trexio_backend), parameter :: TREXIO_HDF5 = 0
integer(trexio_backend), parameter :: TREXIO_TEXT = 1
! integer(trexio_backend), parameter :: TREXIO_JSON = 2
integer(trexio_backend), parameter :: TREXIO_INVALID_BACK_END = 2
#+end_src
#+end_src
@ -584,9 +585,10 @@ trexio_open(const char* file_name, const char mode,
interface
integer(8) function trexio_open_c (filename, mode, backend) bind(C, name="trexio_open")
use, intrinsic :: iso_c_binding
import
character(kind=c_char), dimension(*) :: filename
character, intent(in), value :: mode
integer, intent(in), value :: backend
integer(trexio_backend), intent(in), value :: backend
end function trexio_open_c
end interface
#+end_src
@ -1535,7 +1537,7 @@ contains
implicit none
character(len=*) :: filename
character, intent(in), value :: mode
integer, intent(in), value :: backend
integer(trexio_backend), intent(in), value :: backend
character(len=len_trim(filename)+1) :: filename_c
filename_c = trim(filename) // c_null_char

View File

@ -15,11 +15,15 @@ int main() {
assert (rc == 0);
test_write("test_write.h5", TREXIO_HDF5);
test_read ("test_write.h5", TREXIO_HDF5);
rc = system("rm -rf test_write.h5");
assert (rc == 0);
rc = system("rm -rf test_write.dir");
assert (rc == 0);
test_write("test_write.dir", TREXIO_TEXT);
test_read ("test_write.dir", TREXIO_TEXT);
rc = system("rm -rf test_write.dir");
assert (rc == 0);
return 0;
}

View File

@ -1,21 +1,39 @@
program test_trexio
use trexio
implicit none
call test_write()
call test_read()
call system('rm -rf trexio_test_fort')
print *, 'call test_write(''trexio_test_fort'', TREXIO_TEXT)'
call test_write('trexio_test_fort', TREXIO_TEXT)
print *, 'call test_read(''trexio_test_fort'', TREXIO_TEXT)'
call test_read('trexio_test_fort', TREXIO_TEXT)
call system('rm -rf trexio_test_fort')
call system('rm -rf trexio_test_fort')
print *, 'call test_write(''trexio_test_fort.h5'', TREXIO_HDF5)'
call test_write('trexio_test_fort.h5', TREXIO_HDF5)
print *, 'call test_read(''trexio_test_fort.h5'', TREXIO_HDF5)'
call test_read('trexio_test_fort.h5', TREXIO_HDF5)
call system('rm -rf trexio_test_fort.h5')
end program test_trexio
subroutine test_write()
subroutine test_write(file_name, back_end)
! ============ Test write functionality =============== !
use trexio
implicit none
character*(*), intent(in) :: file_name
integer(trexio_backend), intent(in) :: back_end
integer(8) :: trex_file
integer :: rc = 1
integer :: num
character*(128) :: str
double precision :: charge(12)
double precision :: coord(3,12)
@ -39,33 +57,84 @@ subroutine test_write()
! ================= START OF TEST ===================== !
trex_file = trexio_open('trexio_test_fort', 'w', TREXIO_TEXT)
! trex_file = trexio_open('test_hdf5_fort.h5', 'w', TREXIO_HDF5)
trex_file = trexio_open(file_name, 'w', back_end)
rc = trexio_has_nucleus_num(trex_file)
if (rc == TREXIO_HAS_NOT) write(*,*) 'SUCCESS HAS NOT 1'
if (rc == TREXIO_HAS_NOT) then
write(*,*) 'SUCCESS HAS NOT 1'
else
call trexio_string_of_error(TREXIO_READONLY,str)
print *, trim(str)
call exit(1)
endif
rc = trexio_has_nucleus_charge(trex_file)
if (rc == TREXIO_HAS_NOT) write(*,*) 'SUCCESS HAS NOT 2'
if (rc == TREXIO_HAS_NOT) then
write(*,*) 'SUCCESS HAS NOT 2'
else
call trexio_string_of_error(TREXIO_READONLY,str)
print *, trim(str)
call exit(1)
endif
rc = trexio_write_nucleus_num(trex_file, num)
if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS WRITE NUM'
if (rc == TREXIO_SUCCESS) then
write(*,*) 'SUCCESS WRITE NUM'
else
call trexio_string_of_error(TREXIO_READONLY,str)
print *, trim(str)
call exit(1)
endif
rc = trexio_write_nucleus_charge(trex_file, charge)
if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS WRITE CHARGE'
if (rc == TREXIO_SUCCESS) then
write(*,*) 'SUCCESS WRITE CHARGE'
else
call trexio_string_of_error(TREXIO_READONLY,str)
print *, trim(str)
call exit(1)
endif
rc = trexio_write_nucleus_coord(trex_file, coord)
if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS WRITE COORD'
if (rc == TREXIO_SUCCESS) then
write(*,*) 'SUCCESS WRITE COORD'
else
call trexio_string_of_error(TREXIO_READONLY,str)
print *, trim(str)
call exit(1)
endif
rc = trexio_has_nucleus_num(trex_file)
if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS HAS 1'
if (rc == TREXIO_SUCCESS) then
write(*,*) 'SUCCESS HAS 1'
else
call trexio_string_of_error(TREXIO_READONLY,str)
print *, trim(str)
call exit(1)
endif
rc = trexio_has_nucleus_coord(trex_file)
if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS HAS 2'
if (rc == TREXIO_SUCCESS) then
write(*,*) 'SUCCESS HAS 2'
else
call trexio_string_of_error(TREXIO_READONLY,str)
print *, trim(str)
call exit(1)
endif
rc = trexio_close(trex_file)
if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS CLOSE'
if (rc == TREXIO_SUCCESS) then
write(*,*) 'SUCCESS CLOSE'
else
call trexio_string_of_error(TREXIO_READONLY,str)
print *, trim(str)
call exit(1)
endif
! ---------------------------------- !
! to modify fiels of existing file:
! to modify fields of existing file:
! text backend -> open with 'w'
! hdf5 backend -> open with 'a'
! ---------------------------------- !
@ -85,13 +154,19 @@ subroutine test_write()
end subroutine test_write
subroutine test_read()
subroutine test_read(file_name, back_end)
! ============ Test read functionality =============== !
use trexio
implicit none
character*(*), intent(in) :: file_name
integer(trexio_backend), intent(in) :: back_end
integer(8) :: trex_file
integer :: rc = 1
@ -106,26 +181,48 @@ subroutine test_read()
! ================= START OF TEST ===================== !
trex_file = trexio_open('trexio_test_fort', 'r', TREXIO_TEXT)
! trex_file = trexio_open('test_hdf5_fort.h5', 'r', TREXIO_HDF5)
trex_file = trexio_open(file_name, 'r', back_end)
rc = trexio_read_nucleus_num(trex_file, num_read)
if (rc == TREXIO_SUCCESS .and. num_read == num) write(*,*) 'SUCCESS READ NUM'
if (rc == TREXIO_SUCCESS .and. num_read == num) then
write(*,*) 'SUCCESS READ NUM'
else
call trexio_string_of_error(TREXIO_READONLY,str)
print *, trim(str)
call exit(1)
endif
rc = trexio_read_nucleus_charge(trex_file, charge)
if (rc == TREXIO_SUCCESS .and. (abs(charge(11) - 1.0) < 1.0D-8) ) write(*,*) 'SUCCESS READ CHARGE'
if (rc == TREXIO_SUCCESS .and. (dabs(charge(11) - 1.d0) < 1.0D-8) ) then
write(*,*) 'SUCCESS READ CHARGE'
else
call trexio_string_of_error(TREXIO_READONLY,str)
print *, trim(str)
call exit(-1)
endif
rc = trexio_read_nucleus_coord(trex_file, coord)
if (rc == TREXIO_SUCCESS .and. (abs(coord(2,1) - 1.39250319d0) < 1.0D-8) ) write(*,*) 'SUCCESS READ COORD'
if (rc == TREXIO_SUCCESS .and. (dabs(coord(2,1) - 1.39250319d0) < 1.0D-8) ) then
write(*,*) 'SUCCESS READ COORD'
else
call trexio_string_of_error(TREXIO_READONLY,str)
print *, trim(str)
call exit(-1)
endif
rc = trexio_close(trex_file)
if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS CLOSE'
if (rc == TREXIO_SUCCESS) then
write(*,*) 'SUCCESS CLOSE'
else
call trexio_string_of_error(TREXIO_READONLY,str)
write(*,*) str
print *, trim(str)
call exit(1)
endif
! ================= END OF TEST ===================== !