diff options
author | Elena Pourmal <epourmal@hdfgroup.org> | 2005-09-07 22:24:18 (GMT) |
---|---|---|
committer | Elena Pourmal <epourmal@hdfgroup.org> | 2005-09-07 22:24:18 (GMT) |
commit | 64fdc49ab940e76e921f4ce7ccddc577fb24bea6 (patch) | |
tree | 24f6495fc8cbb1058f1b4f5a8630b41345069a30 /fortran/src | |
parent | ae138025bb35212b021bf0c95f89fb5c51b56448 (diff) | |
download | hdf5-64fdc49ab940e76e921f4ce7ccddc577fb24bea6.zip hdf5-64fdc49ab940e76e921f4ce7ccddc577fb24bea6.tar.gz hdf5-64fdc49ab940e76e921f4ce7ccddc577fb24bea6.tar.bz2 |
[svn-r11369] Purpose: Improvement/maintenance
Description: Added code to generate sizes of Fortran REAL and
DOUBLE PRECISION types. This will "almost" eliminate
H5f90i.h file that defines C stubs datatypes.
Solution:
Platforms tested: heping with g95 (-r8, -d8 and default settings)
Misc. update:
Diffstat (limited to 'fortran/src')
-rw-r--r-- | fortran/src/H5test_kind.f90 | 43 |
1 files changed, 43 insertions, 0 deletions
diff --git a/fortran/src/H5test_kind.f90 b/fortran/src/H5test_kind.f90 index 1f5dca3..7cdc953 100644 --- a/fortran/src/H5test_kind.f90 +++ b/fortran/src/H5test_kind.f90 @@ -1,9 +1,11 @@ ! H5test_kind.f90 ! ! This fortran program generates H5fortran_detect.f90 +! ! program test_kind integer :: i, j, ii, last, kind_numbers(10) + integer :: jr, jd last = -1 ii = 0 j = selected_int_kind(18) @@ -25,6 +27,10 @@ write(*,*) "write(*,*) "" /*generating header file*/ """ j = 0 write(*, "("" call i"", i2.2,""()"")") j + jr = 0 + write(*, "("" call r"", i2.2,""()"")") jr + jd = 0 + write(*, "("" call d"", i2.2,""()"")") jd do i = 1, ii j = kind_numbers(i) write(*, "("" call i"", i2.2,""()"")") j @@ -53,6 +59,43 @@ write(*,*)" endif" write(*,*)" return" write(*,*)" end subroutine" + jr = 0 + write(*, "("" subroutine r"" i2.2,""()"")") j + write(*,*)" implicit none" + write(*,*)" real :: b" + write(*,*)" integer :: a(8)" + write(*,*)" integer :: a_size" + write(*,*)" integer :: b_size" + write(*,*)" a_size = bit_size(a(1))" + write(*,*)" b_size = size(transfer(b,a))*a_size" + write(*,*)" if (b_size .eq. 32) then" + write(*,*)" write(*,*) ""#define H5_FORTRAN_HAS_REAL_NATIVE_4"" " + write(*,*)" endif" + write(*,*)" if (b_size .eq. 64) then" + write(*,*)" write(*,*) ""#define H5_FORTRAN_HAS_REAL_NATIVE_8"" " + write(*,*)" endif" + write(*,*)" if (b_size .eq. 128) then" + write(*,*)" write(*,*) ""#define H5_FORTRAN_HAS_REAL_NATIVE_16"" " + write(*,*)" endif" + write(*,*)" return" + write(*,*)" end subroutine" + jd = 0 + write(*, "("" subroutine d"" i2.2,""()"")") jd + write(*,*)" implicit none" + write(*,*)" double precision :: b" + write(*,*)" integer :: a(8)" + write(*,*)" integer :: a_size" + write(*,*)" integer :: b_size" + write(*,*)" a_size = bit_size(a(1))" + write(*,*)" b_size = size(transfer(b,a))*a_size" + write(*,*)" if (b_size .eq. 64) then" + write(*,*)" write(*,*) ""#define H5_FORTRAN_HAS_DOUBLE_NATIVE_8"" " + write(*,*)" endif" + write(*,*)" if (b_size .eq. 128) then" + write(*,*)" write(*,*) ""#define H5_FORTRAN_HAS_DOUBLE_NATIVE_16"" " + write(*,*)" endif" + write(*,*)" return" + write(*,*)" end subroutine" do i = 1, ii j = kind_numbers(i) write(*, "("" subroutine i"" i2.2,""()"")") j |