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

more fortran-ic trexio_open

This commit is contained in:
q-posev 2021-03-24 17:13:53 +01:00
parent c0b6f95e52
commit f8a28a037d
2 changed files with 34 additions and 16 deletions

View File

@ -295,12 +295,12 @@ trexio_t* trexio_open(const char* file_name, const char mode, const back_end_t b
#+begin_src fortran :tangle prefix_fortran.f90
interface
integer(8) function trexio_open (filename, mode, backend) bind(C)
integer(8) function trexio_open_c (filename, mode, backend) bind(C, name="trexio_open")
use, intrinsic :: iso_c_binding
character(kind=c_char), dimension(*) :: filename
character, intent(in), value :: mode
integer, intent(in), value :: backend
end function trexio_open
character(kind=c_char), dimension(*) :: filename
character, intent(in), value :: mode
integer, intent(in), value :: backend
end function trexio_open_c
end interface
#+end_src
@ -621,6 +621,25 @@ end interface
- Text files: not to be used for production, but useful for debugging
- JSON: for portability
* Fortran helper/wrapper functions
#+begin_src fortran :tangle suffix_fortran.f90
contains
integer(8) function trexio_open (filename, mode, backend)
use, intrinsic :: iso_c_binding
implicit none
character(len=*) :: filename
character, intent(in), value :: mode
integer, intent(in), value :: backend
character(len=len_trim(filename)+1) :: filename_c
filename_c = trim(filename) // c_null_char
trexio_open = trexio_open_c(filename_c, mode, backend)
end function trexio_open
#+end_src
* File suffixes :noxport:
#+begin_src c :tangle suffix_front.h

View File

@ -38,9 +38,8 @@ subroutine test_write()
2.14171677 , 1.23652075 , 0.00000000 , &
0.00000000 , 2.47304151 , 0.00000000 /)
trex_file = trexio_open('test_text_fort' // c_null_char, 'w', TREXIO_TEXT)
! trex_file = trexio_open('test_hdf5_fort.h5' // c_null_char, 'w', TREXIO_HDF5)
! trex_file = trexio_open('test_text_fort', 'w', TREXIO_TEXT)
trex_file = trexio_open('test_hdf5_fort.h5', 'w', TREXIO_HDF5)
rc = trexio_write_nucleus_num(trex_file, num)
if (rc == 0) write(*,*) 'SUCCESS WRITE NUM'
@ -60,16 +59,16 @@ subroutine test_write()
! hdf5 backend -> open with 'a'
! ---------------------------------- !
trex_file = trexio_open('test_text_fort' // c_null_char, 'w', TREXIO_TEXT);
! trex_file = trexio_open('test_text_fort' // c_null_char, 'w', TREXIO_TEXT);
! trex_file = trexio_open('test_hdf5_fort.h5' // c_null_char, 'a', TREXIO_HDF5)
coord(1) = 666.666
! coord(1) = 666.666
rc = trexio_write_nucleus_coord(trex_file,coord)
if (rc == 0) write(*,*) 'SUCCESS MODIFY COORD'
! rc = trexio_write_nucleus_coord(trex_file,coord)
! if (rc == 0) write(*,*) 'SUCCESS MODIFY COORD'
rc = trexio_close(trex_file)
if (rc == 0) write(*,*) 'SUCCESS CLOSE'
! rc = trexio_close(trex_file)
! if (rc == 0) write(*,*) 'SUCCESS CLOSE'
end subroutine test_write
@ -90,8 +89,8 @@ subroutine test_read()
rc = 0
trex_file = trexio_open('test_text_fort' // c_null_char, 'r', TREXIO_TEXT)
! trex_file = trexio_open('test_hdf5_fort.h5' // c_null_char, 'r', TREXIO_HDF5)
! trex_file = trexio_open('test_text_fort', 'r', TREXIO_TEXT)
trex_file = trexio_open('test_hdf5_fort.h5', 'r', TREXIO_HDF5)
rc = trexio_read_nucleus_num(trex_file, num_read)