summaryrefslogtreecommitdiffstats
path: root/fortran/src/H5test_kind.f90
blob: 7cdc95398e0332646d6c2ca0967f6ec8450635b9 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
! 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)
!         write(*,*) j
         do i = 1,100
            j = selected_int_kind(i)
            if(j .ne. last) then
              if(last .ne. -1) then
                  ii = ii + 1
                  kind_numbers(ii) = last
              endif
            last = j
            if(j .eq. -1) exit
            endif
          enddo
!          write(*,*) kind_numbers(1:ii)
! Generate a program
          write(*,*) "program int_kind"
          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
          enddo
          write(*,*) "end program int_kind"
              j = 0
             write(*, "("" subroutine i"" i2.2,""()"")") j
             write(*,*)"   implicit none"
             write(*,*)"   integer :: a"
             write(*,*)"   integer :: a_size"
             write(*,*)"   a_size = bit_size(a)"
             write(*,*)"   if (a_size .eq. 8) then"
             write(*,*)"       write(*,*) ""#define H5_FORTRAN_HAS_NATIVE_1"" "
             write(*,*)"   endif"
             write(*,*)"   if (a_size .eq. 16) then"
             write(*,*)"       write(*,*) ""#define H5_FORTRAN_HAS_NATIVE_2"" "
             write(*,*)"   endif"
             write(*,*)"   if (a_size .eq. 32) then"
             write(*,*)"       write(*,*) ""#define H5_FORTRAN_HAS_NATIVE_4"" "
             write(*,*)"   endif"
             write(*,*)"   if (a_size .eq. 64) then"
             write(*,*)"       write(*,*) ""#define H5_FORTRAN_HAS_NATIVE_8"" "
             write(*,*)"   endif"
             write(*,*)"   if (a_size .eq. 128) then"
             write(*,*)"       write(*,*) ""#define H5_FORTRAN_HAS_NATIVE_16"" "
             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
             write(*,*)"   implicit none"
             write(*,*)"   integer(",j,") :: a"
             write(*,*)"   integer :: a_size"
             write(*,*)"   a_size = bit_size(a)"
             write(*,*)"   if (a_size .eq. 8) then"
             write(*,*)"       write(*,*) ""#define H5_FORTRAN_HAS_INTEGER_1"" "
             write(*,*)"   endif"
             write(*,*)"   if (a_size .eq. 16) then"
             write(*,*)"       write(*,*) ""#define H5_FORTRAN_HAS_INTEGER_2"" "
             write(*,*)"   endif"
             write(*,*)"   if (a_size .eq. 32) then"
             write(*,*)"       write(*,*) ""#define H5_FORTRAN_HAS_INTEGER_4"" "
             write(*,*)"   endif"
             write(*,*)"   if (a_size .eq. 64) then"
             write(*,*)"       write(*,*) ""#define H5_FORTRAN_HAS_INTEGER_8"" "
             write(*,*)"   endif"
             write(*,*)"   if (a_size .eq. 128) then"
             write(*,*)"       write(*,*) ""#define H5_FORTRAN_HAS_INTEGER_16"" "
             write(*,*)"   endif"
             write(*,*)"   return"
             write(*,*)"   end subroutine"    
          enddo
          end program