1
0
mirror of https://github.com/TREX-CoE/trexio.git synced 2024-11-03 20:54:07 +01:00

fix tests

This commit is contained in:
q-posev 2021-06-04 18:42:52 +02:00
parent 1ce1872916
commit 8fb44d585c
2 changed files with 9 additions and 26 deletions

View File

@ -92,7 +92,7 @@ int test_write(const char* file_name, const back_end_t backend) {
if (backend == TREXIO_HDF5) rc = trexio_write_nucleus_label(file,labelxxx, 4);
assert (rc == TREXIO_SUCCESS);
if (backend == TREXIO_HDF5) rc = trexio_write_nucleus_symmetry(file, sym);
if (backend == TREXIO_HDF5) rc = trexio_write_nucleus_point_group(file, sym);
assert (rc == TREXIO_SUCCESS);
// check if the written data exists in the file
rc = trexio_has_nucleus_num(file);
@ -101,7 +101,7 @@ int test_write(const char* file_name, const back_end_t backend) {
assert (rc == TREXIO_SUCCESS);
if (backend == TREXIO_HDF5) rc = trexio_has_nucleus_label(file);
assert (rc == TREXIO_SUCCESS);
if (backend == TREXIO_HDF5) rc = trexio_has_nucleus_symmetry(file);
if (backend == TREXIO_HDF5) rc = trexio_has_nucleus_point_group(file);
assert (rc == TREXIO_SUCCESS);
// should not work: try to overwrite the num variable
@ -146,7 +146,7 @@ int test_read(const char* file_name, const back_end_t backend) {
double* coord;
char** label;
char* labelxxx;
char* symmetry;
char* point_group;
/*================= START OF TEST ==================*/
@ -189,13 +189,13 @@ int test_read(const char* file_name, const back_end_t backend) {
pch = strtok(NULL, TREXIO_DELIM);
assert( strcmp(pch, "Na") == 0 );
symmetry = (char*) malloc(32*sizeof(char));
point_group = (char*) malloc(32*sizeof(char));
rc = trexio_read_nucleus_symmetry(file, symmetry);
rc = trexio_read_nucleus_point_group(file, point_group);
assert (rc == TREXIO_SUCCESS);
assert( strcmp(symmetry, "B3U") == 0 );
free(symmetry);
assert( strcmp(point_group, "B3U") == 0 );
free(point_group);
for (int i=0; i<num; i++){
free(label[i]);

View File

@ -74,7 +74,7 @@ subroutine test_write()
if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS WRITE LABEL'
deallocate(label_str)
rc = trexio_write_nucleus_symmetry(trex_file, sym_str)
rc = trexio_write_nucleus_point_group(trex_file, sym_str)
if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS WRITE SYMMETRY'
deallocate(sym_str)
@ -86,23 +86,6 @@ subroutine test_write()
rc = trexio_close(trex_file)
if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS CLOSE'
! ---------------------------------- !
! to modify fiels of existing file:
! text backend -> open with 'w'
! hdf5 backend -> open with 'a'
! ---------------------------------- !
!! trex_file = trexio_open('trexio_test_fort', 'w', TREXIO_TEXT);
!! trex_file = trexio_open('test_hdf5_fort.h5', 'a', TREXIO_HDF5)
! coord(1) = 666.666
! rc = trexio_write_nucleus_coord(trex_file,coord)
! if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS MODIFY COORD'
! rc = trexio_close(trex_file)
! if (rc == TREXIO_SUCCESS) write(*,*) 'SUCCESS CLOSE'
! ================= END OF TEST ===================== !
end subroutine test_write
@ -177,7 +160,7 @@ subroutine test_read()
if (rc == TREXIO_SUCCESS .and. (trim(label(2)) == 'Na') ) write(*,*) 'SUCCESS READ LABEL'
rc = trexio_read_nucleus_symmetry(trex_file, sym_str)
rc = trexio_read_nucleus_point_group(trex_file, sym_str)
write(*,*) sym_str
if (rc == TREXIO_SUCCESS .and. (trim(sym_str) == 'B3U') ) write(*,*) 'SUCCESS READ SYMMETRY'