1
0
mirror of https://github.com/TREX-CoE/trexio.git synced 2024-08-25 06:31:43 +02:00

Merge pull request #77 from TREX-CoE/iso_fortran_env

Using ISO_C_BINDING types in trexio_f.f90
This commit is contained in:
Evgeny Posenitskiy 2022-01-21 12:32:14 +01:00 committed by GitHub
commit ba2e0691ba
No known key found for this signature in database
GPG Key ID: 4AEE18F83AFDEB23
3 changed files with 151 additions and 105 deletions

View File

@ -25,9 +25,9 @@ cat prefix_python.py > trexio.py
# append version string and attributes to the Fortran module file
echo "" >> trexio_f.f90
echo "character(len = 12) :: TREXIO_PACKAGE_VERSION = ${VERSION_VAL}" >> trexio_f.f90
echo "integer(4) :: TREXIO_VERSION_MAJOR = ${VERSION_MAJOR_VAL}" >> trexio_f.f90
echo "integer(4) :: TREXIO_VERSION_MINOR = ${VERSION_MINOR_VAL}" >> trexio_f.f90
echo "integer(4) :: TREXIO_VERSION_PATCH = ${VERSION_PATCH_VAL}" >> trexio_f.f90
echo "integer :: TREXIO_VERSION_MAJOR = ${VERSION_MAJOR_VAL}" >> trexio_f.f90
echo "integer :: TREXIO_VERSION_MINOR = ${VERSION_MINOR_VAL}" >> trexio_f.f90
echo "integer :: TREXIO_VERSION_PATCH = ${VERSION_PATCH_VAL}" >> trexio_f.f90
echo "" >> trexio_f.f90
# c front end
@ -48,4 +48,3 @@ cat populated/pop_*.py >> trexio.py
cat suffix_s_front.h >> trexio_s.h
cat suffix_front.h >> trexio.h
cat suffix_fortran.f90 >> trexio_f.f90

View File

@ -82,8 +82,9 @@ module trexio
use, intrinsic :: iso_c_binding
implicit none
integer, parameter :: trexio_exit_code = c_int32_t
integer, parameter :: trexio_backend = c_int32_t
integer, parameter :: trexio_exit_code = c_int32_t
integer, parameter :: trexio_back_end_t = c_int32_t
integer, parameter :: trexio_t = c_int64_t
character(kind=c_char), parameter :: TREXIO_DELIM = c_new_line
#+end_src
@ -559,15 +560,16 @@ typedef int32_t back_end_t;
#define TREXIO_DELIM "\n"
#+end_src
The helper function ~trexio_has_backend~ returns ~true~ if TREXIO compilation includes a back end provided as an argument; ~false~ otherwise.
The helper function ~trexio_has_back_end~ returns ~true~ if TREXIO compilation includes a back end provided as an argument; ~false~ otherwise.
This is useful due to the fact that HDF5 back end can be disabled at configure step.
#+begin_src c :tangle prefix_front.h
bool trexio_has_backend(back_end_t back_end);
bool trexio_has_back_end(back_end_t back_end);
#+end_src
#+begin_src c :tangle prefix_front.c
bool trexio_has_backend(back_end_t back_end) {
bool trexio_has_back_end(back_end_t back_end) {
switch (back_end) {
case TREXIO_TEXT:
return true;
@ -580,28 +582,44 @@ bool trexio_has_backend(back_end_t back_end) {
}
return false;
}
bool trexio_has_backend(back_end_t back_end) {
return trexio_has_back_end(back_end);
}
#+end_src
*** Fortran
#+begin_src f90 :tangle prefix_fortran.f90
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
integer(trexio_back_end_t), parameter :: TREXIO_HDF5 = 0
integer(trexio_back_end_t), parameter :: TREXIO_TEXT = 1
! integer(trexio_back_end_t), parameter :: TREXIO_JSON = 2
integer(trexio_back_end_t), parameter :: TREXIO_INVALID_BACK_END = 2
#+end_src
The function below is a Fortran interface for the aforementioned C-compatible ~trexio_has_backend~ function.
The function below is a Fortran interface for the aforementioned C-compatible ~trexio_has_back_end~ function.
#+begin_src f90 :tangle prefix_fortran.f90
interface
logical(c_bool) function trexio_has_back_end (back_end) bind(C)
use, intrinsic :: iso_c_binding
import
integer(trexio_back_end_t), intent(in), value :: back_end
end function trexio_has_back_end
end interface
interface
logical(c_bool) function trexio_has_backend (back_end) bind(C)
use, intrinsic :: iso_c_binding
integer(c_int32_t), intent(in), value :: back_end
import
integer(trexio_back_end_t), intent(in), value :: back_end
end function trexio_has_backend
end interface
#+end_src
Originally, the function was named ~trexio_has_backend~. For
consistency, in version 2.2 it was renamed ~trexio_has_back_end~.
*** Python
#+begin_src python :tangle prefix_python.py
@ -966,13 +984,13 @@ trexio_open(const char* file_name, const char mode,
#+begin_src f90 :tangle prefix_fortran.f90
interface
integer(c_int64_t) function trexio_open_c (filename, mode, backend, rc_open) bind(C, name="trexio_open")
integer(trexio_t) function trexio_open_c (filename, mode, back_end, rc_open) bind(C, name="trexio_open")
use, intrinsic :: iso_c_binding
import
character(kind=c_char), dimension(*) :: filename
character(kind=c_char), intent(in), value :: mode
integer(trexio_backend), intent(in), value :: backend
integer(trexio_exit_code), intent(out) :: rc_open
character(kind=c_char), dimension(*) :: filename
character(kind=c_char), intent(in), value :: mode
integer(trexio_back_end_t), intent(in), value :: back_end
integer(trexio_exit_code), intent(out) :: rc_open
end function trexio_open_c
end interface
#+end_src
@ -1044,8 +1062,9 @@ trexio_exit_code trexio_set_one_based(trexio_t* file)
#+begin_src f90 :tangle prefix_fortran.f90
interface
integer(c_int32_t) function trexio_set_one_based(trex_file) bind(C)
integer(trexio_exit_code) function trexio_set_one_based(trex_file) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
end function trexio_set_one_based
end interface
@ -1147,8 +1166,9 @@ trexio_close (trexio_t* file)
#+begin_src f90 :tangle prefix_fortran.f90
interface
integer(c_int32_t) function trexio_close (trex_file) bind(C)
integer(trexio_exit_code) function trexio_close (trex_file) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
end function trexio_close
end interface
@ -1219,7 +1239,7 @@ trexio_inquire (const char* file_name)
#+begin_src f90 :tangle prefix_fortran.f90
interface
integer function trexio_inquire_c (filename) bind(C, name="trexio_inquire")
integer(trexio_exit_code) function trexio_inquire_c (filename) bind(C, name="trexio_inquire")
use, intrinsic :: iso_c_binding
import
character(kind=c_char), dimension(*) :: filename
@ -1344,7 +1364,7 @@ def _inquire(file_name: str) -> bool:
handle. Second parameter is the variable to be written/read
to/from the ~TREXIO~ file (except for ~trexio_has_~ functions).
Suffixes ~_32~ and ~_64~ correspond to API calls dealing with
single and double precision, respectively. The basic
single and real(c_double), respectively. The basic
(non-suffixed) API call on dimensioning variables deals with single
precision (see Table above).
@ -1546,8 +1566,9 @@ trexio_has_$group_num$ (trexio_t* const file)
#+begin_src f90 :tangle write_attr_num_64_front_fortran.f90
interface
integer(c_int32_t) function trexio_write_$group_num$_64 (trex_file, num) bind(C)
integer(trexio_exit_code) function trexio_write_$group_num$_64 (trex_file, num) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_num_f_dtype_double$, intent(in), value :: num
end function trexio_write_$group_num$_64
@ -1556,8 +1577,9 @@ end interface
#+begin_src f90 :tangle read_attr_num_64_front_fortran.f90
interface
integer(c_int32_t) function trexio_read_$group_num$_64 (trex_file, num) bind(C)
integer(trexio_exit_code) function trexio_read_$group_num$_64 (trex_file, num) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_num_f_dtype_double$, intent(out) :: num
end function trexio_read_$group_num$_64
@ -1566,8 +1588,9 @@ end interface
#+begin_src f90 :tangle write_attr_num_32_front_fortran.f90
interface
integer(c_int32_t) function trexio_write_$group_num$_32 (trex_file, num) bind(C)
integer(trexio_exit_code) function trexio_write_$group_num$_32 (trex_file, num) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_num_f_dtype_single$, intent(in), value :: num
end function trexio_write_$group_num$_32
@ -1576,8 +1599,9 @@ end interface
#+begin_src f90 :tangle read_attr_num_32_front_fortran.f90
interface
integer(c_int32_t) function trexio_read_$group_num$_32 (trex_file, num) bind(C)
integer(trexio_exit_code) function trexio_read_$group_num$_32 (trex_file, num) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_num_f_dtype_single$, intent(out) :: num
end function trexio_read_$group_num$_32
@ -1586,8 +1610,9 @@ end interface
#+begin_src f90 :tangle write_attr_num_def_front_fortran.f90
interface
integer(c_int32_t) function trexio_write_$group_num$ (trex_file, num) bind(C)
integer(trexio_exit_code) function trexio_write_$group_num$ (trex_file, num) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_num_f_dtype_default$, intent(in), value :: num
end function trexio_write_$group_num$
@ -1596,8 +1621,9 @@ end interface
#+begin_src f90 :tangle read_attr_num_def_front_fortran.f90
interface
integer(c_int32_t) function trexio_read_$group_num$ (trex_file, num) bind(C)
integer(trexio_exit_code) function trexio_read_$group_num$ (trex_file, num) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_num_f_dtype_default$, intent(out) :: num
end function trexio_read_$group_num$
@ -1606,8 +1632,9 @@ end interface
#+begin_src f90 :tangle has_attr_num_front_fortran.f90
interface
integer(c_int32_t) function trexio_has_$group_num$ (trex_file) bind(C)
integer(trexio_exit_code) function trexio_has_$group_num$ (trex_file) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
end function trexio_has_$group_num$
end interface
@ -1715,7 +1742,7 @@ def has_$group_num$(trexio_file) -> bool:
First parameter is the ~TREXIO~ file handle. Second parameter is the variable to be written/read
to/from the ~TREXIO~ file (except for ~trexio_has_~ functions).
Suffixes ~_32~ and ~_64~ correspond to API calls dealing with single and double precision, respectively.
The basic (non-suffixed) API call on datasets deals with double precision (see Table above).
The basic (non-suffixed) API call on datasets deals with real(c_double) (see Table above).
**** Function declarations
@ -1999,7 +2026,7 @@ trexio_write_$group_dset$_32 (trexio_t* const file, const $group_dset_dtype_sing
$group_dset_dtype_double$* $group_dset$_64 = CALLOC(dim_size, $group_dset_dtype_double$);
if ($group_dset$_64 == NULL) return TREXIO_ALLOCATION_FAILED;
/* A type conversion from single precision to double reqired since back end only accepts 64-bit data */
/* A type conversion from single precision to double required since back end only accepts 64-bit data */
if ($is_index$) {
for (uint64_t i=0; i<dim_size; ++i){
$group_dset$_64[i] = ($group_dset_dtype_double$) $group_dset$[i] - ($group_dset_dtype_double$) 1;
@ -2209,8 +2236,9 @@ trexio_has_$group_dset$ (trexio_t* const file)
#+begin_src f90 :tangle write_dset_data_64_front_fortran.f90
interface
integer(c_int32_t) function trexio_write_$group_dset$_64 (trex_file, dset) bind(C)
integer(trexio_exit_code) function trexio_write_$group_dset$_64 (trex_file, dset) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_dset_f_dtype_double$, intent(in) :: dset$group_dset_f_dims$
end function trexio_write_$group_dset$_64
@ -2219,8 +2247,9 @@ end interface
#+begin_src f90 :tangle read_dset_data_64_front_fortran.f90
interface
integer(c_int32_t) function trexio_read_$group_dset$_64 (trex_file, dset) bind(C)
integer(trexio_exit_code) function trexio_read_$group_dset$_64 (trex_file, dset) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_dset_f_dtype_double$, intent(out) :: dset$group_dset_f_dims$
end function trexio_read_$group_dset$_64
@ -2229,8 +2258,9 @@ end interface
#+begin_src f90 :tangle write_dset_data_32_front_fortran.f90
interface
integer(c_int32_t) function trexio_write_$group_dset$_32 (trex_file, dset) bind(C)
integer(trexio_exit_code) function trexio_write_$group_dset$_32 (trex_file, dset) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_dset_f_dtype_single$, intent(in) :: dset$group_dset_f_dims$
end function trexio_write_$group_dset$_32
@ -2239,8 +2269,9 @@ end interface
#+begin_src f90 :tangle read_dset_data_32_front_fortran.f90
interface
integer(c_int32_t) function trexio_read_$group_dset$_32 (trex_file, dset) bind(C)
integer(trexio_exit_code) function trexio_read_$group_dset$_32 (trex_file, dset) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_dset_f_dtype_single$, intent(out) :: dset$group_dset_f_dims$
end function trexio_read_$group_dset$_32
@ -2249,8 +2280,9 @@ end interface
#+begin_src f90 :tangle write_dset_data_def_front_fortran.f90
interface
integer(c_int32_t) function trexio_write_$group_dset$ (trex_file, dset) bind(C)
integer(trexio_exit_code) function trexio_write_$group_dset$ (trex_file, dset) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_dset_f_dtype_default$, intent(in) :: dset$group_dset_f_dims$
end function trexio_write_$group_dset$
@ -2259,8 +2291,9 @@ end interface
#+begin_src f90 :tangle read_dset_data_def_front_fortran.f90
interface
integer(c_int32_t) function trexio_read_$group_dset$ (trex_file, dset) bind(C)
integer(trexio_exit_code) function trexio_read_$group_dset$ (trex_file, dset) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
$group_dset_f_dtype_default$, intent(out) :: dset$group_dset_f_dims$
end function trexio_read_$group_dset$
@ -2269,8 +2302,9 @@ end interface
#+begin_src f90 :tangle has_dset_data_front_fortran.f90
interface
integer(c_int32_t) function trexio_has_$group_dset$ (trex_file) bind(C)
integer(trexio_exit_code) function trexio_has_$group_dset$ (trex_file) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
end function trexio_has_$group_dset$
end interface
@ -2818,10 +2852,11 @@ trexio_has_$group_dset$ (trexio_t* const file)
#+begin_src f90 :tangle write_dset_sparse_front_fortran.f90
interface
integer(c_int32_t) function trexio_write_$group_dset$ (trex_file, &
offset_file, buffer_size, &
index_sparse, value_sparse) bind(C)
integer(trexio_exit_code) function trexio_write_$group_dset$ (trex_file, &
offset_file, buffer_size, &
index_sparse, value_sparse) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
integer(c_int64_t), intent(in), value :: offset_file
integer(c_int64_t), intent(in), value :: buffer_size
@ -2831,11 +2866,12 @@ interface
end interface
interface
integer(c_int32_t) function trexio_write_safe_$group_dset$ (trex_file, &
offset_file, buffer_size, &
index_sparse, index_size, &
value_sparse, value_size) bind(C)
integer(trexio_exit_code) function trexio_write_safe_$group_dset$ (trex_file, &
offset_file, buffer_size, &
index_sparse, index_size, &
value_sparse, value_size) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
integer(c_int64_t), intent(in), value :: offset_file
integer(c_int64_t), intent(in), value :: buffer_size
@ -2849,10 +2885,11 @@ end interface
#+begin_src f90 :tangle read_dset_sparse_front_fortran.f90
interface
integer(c_int32_t) function trexio_read_$group_dset$ (trex_file, &
offset_file, buffer_size, &
index_sparse, value_sparse) bind(C)
integer(trexio_exit_code) function trexio_read_$group_dset$ (trex_file, &
offset_file, buffer_size, &
index_sparse, value_sparse) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
integer(c_int64_t), intent(in), value :: offset_file
integer(c_int64_t), intent(inout) :: buffer_size
@ -2862,11 +2899,12 @@ interface
end interface
interface
integer(c_int32_t) function trexio_read_safe_$group_dset$ (trex_file, &
offset_file, buffer_size, &
index_sparse, index_size, &
value_sparse, value_size) bind(C)
integer(trexio_exit_code) function trexio_read_safe_$group_dset$ (trex_file, &
offset_file, buffer_size, &
index_sparse, index_size, &
value_sparse, value_size) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
integer(c_int64_t), intent(in), value :: offset_file
integer(c_int64_t), intent(inout) :: buffer_size
@ -2880,19 +2918,21 @@ end interface
#+begin_src f90 :tangle read_dset_sparse_size_front_fortran.f90
interface
integer(c_int32_t) function trexio_read_$group_dset$_size (trex_file, &
size_max) bind(C)
integer(trexio_exit_code) function trexio_read_$group_dset$_size (trex_file, &
size_max) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
integer(c_int64_t), intent(out) :: size_max
integer(c_int64_t), intent(out) :: size_max
end function trexio_read_$group_dset$_size
end interface
#+end_src
#+begin_src f90 :tangle has_dset_sparse_front_fortran.f90
interface
integer(c_int32_t) function trexio_has_$group_dset$ (trex_file) bind(C)
integer(trexio_exit_code) function trexio_has_$group_dset$ (trex_file) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
end function trexio_has_$group_dset$
end interface
@ -3389,8 +3429,9 @@ trexio_has_$group_dset$ (trexio_t* const file)
#+begin_src f90 :tangle write_dset_str_front_fortran.f90
interface
integer(c_int32_t) function trexio_write_$group_dset$_low (trex_file, dset, max_str_len) bind(C)
integer(trexio_exit_code) function trexio_write_$group_dset$_low (trex_file, dset, max_str_len) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
character(kind=c_char), intent(in) :: dset(*)
integer(c_int32_t), intent(in), value :: max_str_len
@ -3400,8 +3441,9 @@ end interface
#+begin_src f90 :tangle read_dset_str_front_fortran.f90
interface
integer(c_int32_t) function trexio_read_$group_dset$_low (trex_file, dset, max_str_len) bind(C)
integer(trexio_exit_code) function trexio_read_$group_dset$_low (trex_file, dset, max_str_len) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
character(kind=c_char), intent(out) :: dset(*)
integer(c_int32_t), intent(in), value :: max_str_len
@ -3411,23 +3453,24 @@ end interface
#+begin_src f90 :tangle has_dset_str_front_fortran.f90
interface
integer(c_int32_t) function trexio_has_$group_dset$ (trex_file) bind(C)
integer(trexio_exit_code) function trexio_has_$group_dset$ (trex_file) bind(C)
use, intrinsic :: iso_c_binding
import
integer(c_int64_t), intent(in), value :: trex_file
end function trexio_has_$group_dset$
end interface
#+end_src
#+begin_src f90 :tangle helper_read_dset_str_front_fortran.fh_90
integer function trexio_read_$group_dset$ (trex_file, dset, max_str_len)
integer(trexio_exit_code) function trexio_read_$group_dset$ (trex_file, dset, max_str_len)
implicit none
integer(8), intent(in), value :: trex_file
integer(4), intent(in), value :: max_str_len
integer(c_int64_t), intent(in), value :: trex_file
integer(c_int32_t), intent(in), value :: max_str_len
character(len=*), intent(inout) :: dset(*)
character, allocatable :: str_compiled(:)
integer(8) :: $group_dset_dim$
integer :: rc
integer(c_int64_t) :: $group_dset_dim$
integer(trexio_exit_code) :: rc
rc = trexio_read_$group_dset_dim$_64(trex_file, $group_dset_dim$)
if (rc /= TREXIO_SUCCESS) trexio_read_$group_dset$ = rc
@ -3448,15 +3491,15 @@ end interface
#+end_src
#+begin_src f90 :tangle helper_write_dset_str_front_fortran.fh_90
integer function trexio_write_$group_dset$ (trex_file, dset, max_str_len)
integer(trexio_exit_code) function trexio_write_$group_dset$ (trex_file, dset, max_str_len)
implicit none
integer(8), intent(in), value :: trex_file
integer(4), intent(in), value :: max_str_len
character(len=*), intent(in) :: dset(*)
integer(c_int64_t), intent(in), value :: trex_file
integer(c_int32_t), intent(in), value :: max_str_len
character(len=*), intent(in) :: dset(*)
character(len=:), allocatable :: str_compiled
integer(8) :: $group_dset_dim$
integer :: rc
integer(c_int64_t) :: $group_dset_dim$
integer(trexio_exit_code) :: rc
rc = trexio_read_$group_dset_dim$_64(trex_file, $group_dset_dim$)
if (rc /= TREXIO_SUCCESS) then
@ -3705,10 +3748,11 @@ trexio_has_$group_str$ (trexio_t* const file)
#+begin_src f90 :tangle write_attr_str_front_fortran.f90
interface
integer(c_int32_t) function trexio_write_$group_str$_c (trex_file, str, max_str_len) &
bind(C, name="trexio_write_$group_str$")
integer(trexio_exit_code) function trexio_write_$group_str$_c (trex_file, str, max_str_len) &
bind(C, name="trexio_write_$group_str$")
use, intrinsic :: iso_c_binding
integer(c_int64_t), intent(in), value :: trex_file
import
integer(trexio_t), intent(in), value :: trex_file
character(kind=c_char), intent(in) :: str(*)
integer(c_int32_t), intent(in), value :: max_str_len
end function trexio_write_$group_str$_c
@ -3717,30 +3761,32 @@ end interface
#+begin_src f90 :tangle read_attr_str_front_fortran.f90
interface
integer(c_int32_t) function trexio_read_$group_str$_c (trex_file, str, max_str_len) &
bind(C, name="trexio_read_$group_str$")
integer(trexio_exit_code) function trexio_read_$group_str$_c (trex_file, str, max_str_len) &
bind(C, name="trexio_read_$group_str$")
use, intrinsic :: iso_c_binding
integer(c_int64_t), intent(in), value :: trex_file
import
integer(trexio_t), intent(in), value :: trex_file
character(kind=c_char), intent(out) :: str(*)
integer(c_int32_t), intent(in), value :: max_str_len
integer(c_int32_t), intent(in), value :: max_str_len
end function trexio_read_$group_str$_c
end interface
#+end_src
#+begin_src f90 :tangle has_attr_str_front_fortran.f90
interface
integer(c_int32_t) function trexio_has_$group_str$ (trex_file) bind(C)
integer(trexio_exit_code) function trexio_has_$group_str$ (trex_file) bind(C)
use, intrinsic :: iso_c_binding
integer(c_int64_t), intent(in), value :: trex_file
import
integer(trexio_t), intent(in), value :: trex_file
end function trexio_has_$group_str$
end interface
#+end_src
#+begin_src f90 :tangle helper_read_attr_str_front_fortran.fh_90
integer function trexio_read_$group_str$ (trex_file, str, max_str_len)
integer(trexio_exit_code) function trexio_read_$group_str$ (trex_file, str, max_str_len)
implicit none
integer(8), intent(in), value :: trex_file
integer(4), intent(in), value :: max_str_len
integer(trexio_t), intent(in), value :: trex_file
integer(c_int32_t), intent(in), value :: max_str_len
character, intent(out) :: str(*)
trexio_read_$group_str$ = trexio_read_$group_str$_c(trex_file, str, max_str_len)
@ -3749,12 +3795,12 @@ end interface
#+end_src
#+begin_src f90 :tangle helper_write_attr_str_front_fortran.fh_90
integer function trexio_write_$group_str$ (trex_file, str, max_str_len)
integer(trexio_exit_code) function trexio_write_$group_str$ (trex_file, str, max_str_len)
use, intrinsic :: iso_c_binding, only : c_null_char
implicit none
integer(8), intent(in), value :: trex_file
integer(4), intent(in), value :: max_str_len
character(len=*), intent(in) :: str
integer(trexio_t), intent(in), value :: trex_file
integer(c_int32_t), intent(in), value :: max_str_len
character(len=*), intent(in) :: str
character(len=len_trim(str)+1) :: str_c
@ -3857,18 +3903,18 @@ def has_$group_str$(trexio_file) -> bool:
#+begin_src f90 :tangle helper_fortran.f90
contains
integer(8) function trexio_open (filename, mode, backend, rc_open)
integer(trexio_t) function trexio_open (filename, mode, back_end, rc_open)
use, intrinsic :: iso_c_binding, only : c_null_char
implicit none
character(len=*), intent(in) :: filename
character, intent(in), value :: mode
integer(trexio_backend), intent(in), value :: backend
integer(trexio_exit_code), intent(out) :: rc_open
character(len=len_trim(filename)+1) :: filename_c
character(len=*), intent(in) :: filename
character, intent(in), value :: mode
integer(trexio_back_end_t), intent(in), value :: back_end
integer(trexio_exit_code), intent(out) :: rc_open
character(len=len_trim(filename)+1) :: filename_c
integer(trexio_exit_code) :: rc
filename_c = trim(filename) // c_null_char
trexio_open = trexio_open_c(filename_c, mode, backend, rc_open)
trexio_open = trexio_open_c(filename_c, mode, back_end, rc_open)
if (trexio_open == 0_8 .or. rc_open /= TREXIO_SUCCESS) then
return
endif
@ -3885,7 +3931,7 @@ contains
Note, that Fortran interface calls the main ~TREXIO~ API, which is written in C.
#+begin_src f90 :tangle helper_fortran.f90
integer function trexio_inquire (filename)
integer(trexio_exit_code) function trexio_inquire (filename)
use, intrinsic :: iso_c_binding
implicit none
character(len=*), intent(in) :: filename
@ -3906,11 +3952,11 @@ contains
use, intrinsic :: iso_c_binding, only : c_null_char
implicit none
integer(8), intent(in), value :: max_num_str ! number of elements in strign array
integer(c_int64_t), intent(in), value :: max_num_str ! number of elements in strign array
integer, intent(in), value :: max_len_str ! maximum length of a string in an array
character(len=*), intent(in) :: str_array(*)
character(len=:), allocatable, intent(out) :: str_res
integer(8) :: i
integer(c_int64_t) :: i
str_res = ''
do i = 1, max_num_str
@ -3927,13 +3973,13 @@ contains
subroutine trexio_str2strarray(str_flat, max_num_str, max_len_str, str_array)
implicit none
integer(8), intent(in), value :: max_num_str ! number of elements in strign array
integer, intent(in), value :: max_len_str ! maximum length of a string in an array
character, intent(in) :: str_flat(*)
character(len=*), intent(inout) :: str_array(*)
integer(c_int64_t), intent(in), value :: max_num_str ! number of elements in strign array
integer, intent(in), value :: max_len_str ! maximum length of a string in an array
character(kind=c_char), intent(in) :: str_flat(*)
character(len=*), intent(inout) :: str_array(*)
character(len=max_len_str) :: tmp_str
integer(8) :: len_flat, i, j, k, ind
integer(c_int64_t) :: i, j, k, ind, len_flat
len_flat = (max_len_str+1)*max_num_str + 1

View File

@ -384,7 +384,7 @@ exponent =
coefficient =
[ 0.006068, 0.045308, 0.202822, 0.503903, 0.383421, 1.0, 1.0, 1.0, 1.0, 1.0,
0.006068, 0.045308, 0.202822, 0.503903, 0.383421, 1.0, 1.0, 1.0, 1.0, 1.0 ]
`
prim_factor =
[ 1.0006253235944540e+01, 2.4169531573445120e+00, 7.9610924849766440e-01
3.0734305383061117e-01, 1.2929684417481876e-01, 3.0734305383061117e-01,
@ -401,8 +401,9 @@ prim_factor =
construction of all the angular functions of each shell. We
consider two cases for the angular functions: the real-valued
spherical harmonics, and the polynomials in Cartesian coordinates.
In the case of spherical harmonics, the AOs are ordered as
$0, +1, -1, +2, -2, \dots, +m, -m$ and in the case of polynomials we
In the case of real spherical harmonics, the AOs are ordered as
$0, +1, -1, +2, -2, \dots, +m, -m$ (see [[https://en.wikipedia.org/wiki/Table_of_spherical_harmonics#Real_spherical_harmonics][Wikipedia]]).
In the case of polynomials we
impose the canonical (or alphabetical) ordering), i.e
\begin{eqnarray}