rngs.f90
Go to the documentation of this file.
1!#############################################################################
2!# #
3!# fosite - 3D hydrodynamical simulation program #
4!# module: rngs.f90 #
5!# #
6!# Copyright (C) 2015 Manuel Jung <mjung@astrophysik.uni-kiel.de> #
7!# #
8!# This program is free software; you can redistribute it and/or modify #
9!# it under the terms of the GNU General Public License as published by #
10!# the Free Software Foundation; either version 2 of the License, or (at #
11!# your option) any later version. #
12!# #
13!# This program is distributed in the hope that it will be useful, but #
14!# WITHOUT ANY WARRANTY; without even the implied warranty of #
15!# MERCHANTABILITY OR FITNESS FOR A PARTICULAR PURPOSE, GOOD TITLE or #
16!# NON INFRINGEMENT. See the GNU General Public License for more #
17!# details. #
18!# #
19!# You should have received a copy of the GNU General Public License #
20!# along with this program; if not, write to the Free Software #
21!# Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. #
22!# #
23!#############################################################################
24
25!----------------------------------------------------------------------------!
41!----------------------------------------------------------------------------!
42MODULE rngs
43 IMPLICIT NONE
44 !--------------------------------------------------------------------------!
45 PRIVATE
46#if defined(NECSXAURORA)
47 INTEGER, PARAMETER :: i8 = 8
48#else
49 INTEGER, PARAMETER :: i8 = selected_int_kind(18)
50#endif
51 !--------------------------------------------------------------------------!
52 PUBLIC :: &
53 i8, &
54 kiss64, &
57 !--------------------------------------------------------------------------!
58
59CONTAINS
60
61 FUNCTION kiss64(seed) RESULT(res)
62 IMPLICIT NONE
63 !------------------------------------------------------------------------!
64 INTEGER(KIND=I8) :: res
65 INTEGER(KIND=I8), OPTIONAL :: seed
66 !------------------------------------------------------------------------!
67 INTEGER(KIND=I8), save :: x, y, z, c
68 INTEGER(KIND=I8) :: t
69 data x, y, z, c &
70 / 1234567890987654321_i8, &
71 362436362436362436_i8, &
72 1066149217761810_i8, &
73 123456123456123456_i8 /
74 !------------------------------------------------------------------------!
75 IF(PRESENT(seed)) THEN
76 z = seed
77 res = 0
78 ELSE
79 t = ishft(x, 58) + c
80 IF (ishft(x,-63_i8) .EQ. ishft(t,-63_i8)) THEN
81 c = ishft(x, -6) + ishft(x, -63_i8)
82 ELSE
83 c = ishft(x, -6) + 1 - ishft(x+t, -63_i8)
84 END IF
85 x = t + x
86 y = ieor(y, ishft(y,13_i8))
87 y = ieor(y, ishft(y,-17_i8))
88 y = ieor(y, ishft(y,43_i8))
89
90 z = 6906969069_i8 * z + 1234567
91 res = x + y + z
92 END IF
93 END FUNCTION kiss64
94
95 FUNCTION superkiss64() RESULT(x)
96 IMPLICIT NONE
97 !------------------------------------------------------------------------!
98 INTEGER(KIND=I8) :: x
99 !------------------------------------------------------------------------!
100 INTEGER(KIND=I8), save :: i,q(20632),carry,xcng,xs,indx
101 data i,carry,xcng,xs,indx &
102 / 0_i8, &
103 36243678541_i8, &
104 12367890123456_i8,&
105 521288629546311_i8,&
106 20633_i8 /
107 !------------------------------------------------------------------------!
108! Nec SX internal! compiler error:
109! "f90 fatal: Internal error in optimization phase."
110#if !defined(NECSXAURORA)
111 !Initalize if i=0
112 IF(i.EQ.0) THEN
113 DO i=1,20632
114 !fill Q with Congruential+Xorshift
115 xcng = xcng*6906969069_i8+123
116 xs = ieor(xs,ishft(xs,13))
117 xs = ieor(xs,ishft(xs,-17))
118 xs = ieor(xs,ishft(xs,43))
119 q(i) = xcng + xs
120 END DO
121 END IF
122
123 IF(indx <= 20632) THEN
124 x=q(indx)
125 indx=indx+1
126 ELSE
127 x=refill()
128 END IF
129
130 xcng=xcng*6906969069_i8+123
131 xs=ieor(xs,ishft(xs,13))
132 xs=ieor(xs,ishft(xs,-17))
133 xs=ieor(xs,ishft(xs,43))
134 x=x+xcng+xs
135 CONTAINS
136
137 FUNCTION refill() RESULT(s)
138 IMPLICIT NONE
139 !------------------------------------------------------------------------!
140 INTEGER(KIND=I8) :: i,s,z,h
141 !------------------------------------------------------------------------!
142 DO i=1,20632
143 h=iand(carry,1_i8)
144 z = ishft(ishft(q(i),41),-1) &
145 + ishft(ishft(q(i),39),-1) &
146 + ishft(carry,-1)
147 carry=ishft(q(i),-23)+ishft(q(i),-25)+ishft(z,-63)
148 q(i)=not(ishft(z,1)+h)
149 END DO
150 indx=2
151 s=q(1)
152 RETURN
153 END FUNCTION refill
154#else
155 x = 0_i8
156#endif
157 END FUNCTION superkiss64
158
160 FUNCTION i8tod(x) RESULT(res)
161 IMPLICIT NONE
162 !------------------------------------------------------------------------!
163 INTEGER(KIND=I8) :: x
164 REAL :: res
165#if defined(NECSXAURORA)
166 INTEGER(KIND=8) :: long
167#endif
168 !------------------------------------------------------------------------!
169#if defined(NECSXAURORA)
170 ! No idea what is happeing on the SX.., but HUGE(x)!=HUGE(long)
171 ! The 2048? This is 2.*1024. No idea, why the 1024 is needed. This shouldn't
172 ! be required.
173 res = (real(x) / huge(long)) / 2048. + 0.5
174#else
175 res = 0.5e+0 * (real(x) / real(huge(x))) + 0.5e+0
176#endif
177 END FUNCTION i8tod
178
179
180 FUNCTION dkiss64() RESULT(res)
181 IMPLICIT NONE
182 !------------------------------------------------------------------------!
183 REAL :: res
184 !------------------------------------------------------------------------!
185 res = i8tod(kiss64())
186 END FUNCTION dkiss64
187
188
189END MODULE rngs
Definition: rngs.f90:42
real function, public dkiss64()
Definition: rngs.f90:181
integer(kind=i8) function, public kiss64(seed)
Definition: rngs.f90:62
integer(kind=i8) function, public superkiss64()
Definition: rngs.f90:96
integer, parameter, public i8
Definition: rngs.f90:47
real function i8tod(x)
Convert integer I8 random numbers to double precision real values in [0,1].
Definition: rngs.f90:161
integer(kind=i8) function refill()
Definition: rngs.f90:138