A Modular Arbitrary-Order Ocean-Atmosphere Model -- Stochastic implementation
util.f90
Go to the documentation of this file.
1 
2 ! util.f90
3 !
4 !> Utility module
5 !
6 !> @copyright
7 !> 2015 Lesley De Cruz & Jonathan Demaeyer.
8 !> See LICENSE.txt for license information.
9 !
10 !---------------------------------------------------------------------------!
11 
12 MODULE util
13  IMPLICIT NONE
14 
15  PRIVATE
16 
19  PUBLIC :: vector_outer
20 
21 CONTAINS
22 
23  ! SUBROUTINE scalar_allocate(x)
24  ! INTEGER :: AllocStat
25  ! IF (.NOT. ALLOCATED(x)) THEN
26  ! ALLOCATE(x, STAT=AllocStat)
27 
28  !> Convert an integer to string.
29  CHARACTER(len=20) FUNCTION str(k)
30  INTEGER, INTENT(IN) :: k
31  WRITE (str, *) k
32  str = adjustl(str)
33  END FUNCTION str
34 
35  !> Convert a real to string with a given format
36  CHARACTER(len=40) FUNCTION rstr(x,fm)
37  REAL(KIND=8), INTENT(IN) :: x
38  CHARACTER(len=20), INTENT(IN) :: fm
39  WRITE (rstr, trim(adjustl(fm))) x
40  rstr = adjustl(rstr)
41  END FUNCTION rstr
42 
43  !> Random generator initialization routine
44  SUBROUTINE init_random_seed()
45  USE iso_fortran_env, only: int64
46  USE ifport !, only: getpid
47  IMPLICIT NONE
48  INTEGER, ALLOCATABLE :: seed_loc(:)
49  INTEGER :: i, n, un, istat, dt(8), pid
50  INTEGER(int64) :: t
51 
52  CALL random_seed(size = n)
53  ALLOCATE(seed_loc(n))
54  ! First try IF the OS provides a random number generator
55  OPEN(newunit=un, file="/dev/urandom", access="stream", &
56  form="unformatted", action="read", status="old", iostat=istat)
57  IF (istat == 0) THEN
58  READ(un) seed_loc
59  CLOSE(un)
60  ELSE
61  ! Fallback to XOR:ing the current time and pid. The PID is
62  ! useful in case one launches multiple instances of the same
63  ! program in parallel.
64  CALL system_clock(t)
65  IF (t == 0) THEN
66  CALL date_and_time(values=dt)
67  t = (dt(1) - 1970) * 365_int64 * 24 * 60 * 60 * 1000 &
68  + dt(2) * 31_int64 * 24 * 60 * 60 * 1000 &
69  + dt(3) * 24_int64 * 60 * 60 * 1000 &
70  + dt(5) * 60 * 60 * 1000 &
71  + dt(6) * 60 * 1000 + dt(7) * 1000 &
72  + dt(8)
73  END IF
74  pid = getpid()
75  t = ieor(t, int(pid, kind(t)))
76  DO i = 1, n
77  seed_loc(i) = lcg(t)
78  END DO
79  END IF
80  CALL random_seed(put=seed_loc)
81  contains
82  ! This simple PRNG might not be good enough for real work, but is
83  ! sufficient for seeding a better PRNG.
84  FUNCTION lcg(s)
85  integer :: lcg
86  integer(int64) :: s
87  IF (s == 0) THEN
88  s = 104729
89  ELSE
90  s = mod(s, 4294967296_int64)
91  END IF
92  s = mod(s * 279470273_int64, 4294967291_int64)
93  lcg = int(mod(s, int(huge(0), int64)), kind(0))
94  END FUNCTION lcg
95  END SUBROUTINE init_random_seed
96 
97 
98  !> Initialize a square matrix A as a unit matrix
99  SUBROUTINE init_one(A)
100  REAL(KIND=8), DIMENSION(:,:),INTENT(INOUT) :: A
101  INTEGER :: i,n
102  n=size(a,1)
103  a=0.0d0
104  DO i=1,n
105  a(i,i)=1.0d0
106  END DO
107 
108  END SUBROUTINE init_one
109 
110  FUNCTION mat_trace(A)
111  REAL(KIND=8), DIMENSION(:,:) :: A
112  REAL(KIND=8) :: mat_trace
113  INTEGER :: i,n
114  n=size(a,1)
115  mat_trace=0.d0
116  DO i=1,n
117  mat_trace=mat_trace+a(i,i)
118  END DO
119  RETURN
120  END FUNCTION mat_trace
121 
122  FUNCTION mat_contract(A,B)
123  REAL(KIND=8), DIMENSION(:,:) :: A,B
124  REAL(KIND=8) :: mat_contract
125  INTEGER :: i,j,n
126  n=size(a,1)
127  mat_contract=0.d0
128  DO i=1,n
129  DO j=1,n
130  mat_contract=mat_contract+a(i,j)*b(i,j)
131  END DO
132  ENDDO
133  RETURN
134  END FUNCTION mat_contract
135 
136  SUBROUTINE choldc(a,p)
137  REAL(KIND=8), DIMENSION(:,:) :: a
138  REAL(KIND=8), DIMENSION(:) :: p
139  INTEGER :: n
140  INTEGER :: i,j,k
141  REAL(KIND=8) :: sum
142  n=size(a,1)
143  DO i=1,n
144  DO j=i,n
145  sum=a(i,j)
146  DO k=i-1,1,-1
147  sum=sum-a(i,k)*a(j,k)
148  END DO
149  IF (i.eq.j) THEN
150  IF (sum.le.0.) stop 'choldc failed'
151  p(i)=sqrt(sum)
152  ELSE
153  a(j,i)=sum/p(i)
154  ENDIF
155  END DO
156  END DO
157  RETURN
158  END SUBROUTINE choldc
159 
160  SUBROUTINE printmat(A) ! to be moved to util
161  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: A
162  INTEGER :: i
163 
164  DO i=1,SIZE(a,1)
165  print*, a(i,:)
166  END DO
167  END SUBROUTINE printmat
168 
169  SUBROUTINE cprintmat(A) ! to be moved to util
170  COMPLEX(KIND=16), DIMENSION(:,:), INTENT(IN) :: A
171  INTEGER :: i
172 
173  DO i=1,SIZE(a,1)
174  print*, a(i,:)
175  END DO
176  END SUBROUTINE cprintmat
177 
178  FUNCTION invmat(A) RESULT(Ainv)
179  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: A
180  REAL(KIND=8), DIMENSION(SIZE(A,1),SIZE(A,2)) :: Ainv
181 
182  REAL(KIND=8), DIMENSION(SIZE(A,1)) :: work ! work array for LAPACK
183  INTEGER, DIMENSION(SIZE(A,1)) :: ipiv ! pivot indices
184  INTEGER :: n, info
185 
186  ! Store A in Ainv to prevent it from being overwritten by LAPACK
187  ainv = a
188  n = size(a,1)
189 
190  ! DGETRF computes an LU factorization of a general M-by-N matrix A
191  ! using partial pivoting with row interchanges.
192  CALL dgetrf(n, n, ainv, n, ipiv, info)
193 
194  IF (info /= 0) THEN
195  stop 'Matrix is numerically singular!'
196  ENDIF
197 
198  ! DGETRI computes the inverse of a matrix using the LU factorization
199  ! computed by DGETRF.
200  CALL dgetri(n, ainv, n, ipiv, work, n, info)
201 
202  IF (info /= 0) THEN
203  stop 'Matrix inversion failed!'
204  ENDIF
205  END FUNCTION invmat
206 
207  SUBROUTINE triu(A,T)
208  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: A
209  REAL(KIND=8), DIMENSION(:,:), INTENT(OUT) :: T
210  INTEGER i,j
211  t=0.d0
212  DO i=1,SIZE(a,1)
213  DO j=i,SIZE(a,1)
214  t(i,j)=a(i,j)
215  END DO
216  END DO
217  END SUBROUTINE triu
218 
219  SUBROUTINE diag(A,d)
220  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: A
221  REAL(KIND=8), DIMENSION(:), INTENT(OUT) :: d
222  INTEGER :: i
223 
224  DO i=1,SIZE(a,1)
225  d(i)=a(i,i)
226  END DO
227  END SUBROUTINE diag
228 
229  SUBROUTINE cdiag(A,d)
230  COMPLEX(KIND=16), DIMENSION(:,:), INTENT(IN) :: A
231  COMPLEX(KIND=16), DIMENSION(:), INTENT(OUT) :: d
232  INTEGER :: i
233 
234  DO i=1,SIZE(a,1)
235  d(i)=a(i,i)
236  END DO
237  END SUBROUTINE cdiag
238 
239 
240  FUNCTION floordiv(i,j)
241  INTEGER :: i,j,floordiv
242  floordiv=int(floor(real(i)/real(j)))
243  RETURN
244  END FUNCTION floordiv
245 
246  SUBROUTINE reduce(A,Ared,n,ind,rind)
247  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: A
248  REAL(KIND=8), DIMENSION(:,:), INTENT(OUT) :: Ared
249  INTEGER, INTENT(OUT) :: n
250  INTEGER, DIMENSION(:), INTENT(OUT) :: ind,rind
251  LOGICAL, DIMENSION(SIZE(A,1)) :: sel
252  INTEGER :: i,j
253 
254  ind=0
255  rind=0
256  sel=.false.
257  n=0
258  DO i=1,SIZE(a,1)
259  IF (any(a(i,:)/=0)) THEN
260  n=n+1
261  sel(i)=.true.
262  ind(n)=i
263  rind(i)=n
264  ENDIF
265  END DO
266  ared=0.d0
267  DO i=1,SIZE(a,1)
268  DO j=1,SIZE(a,1)
269  IF (sel(i).and.sel(j)) ared(rind(i),rind(j))=a(i,j)
270  ENDDO
271  ENDDO
272  END SUBROUTINE reduce
273 
274  SUBROUTINE ireduce(A,Ared,n,ind,rind)
275  REAL(KIND=8), DIMENSION(:,:), INTENT(OUT) :: A
276  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: Ared
277  INTEGER, INTENT(IN) :: n
278  INTEGER, DIMENSION(:), INTENT(IN) :: ind,rind
279  INTEGER :: i,j
280  a=0.d0
281  DO i=1,n
282  DO j=1,n
283  a(ind(i),ind(j))=ared(i,j)
284  END DO
285  END DO
286  END SUBROUTINE ireduce
287 
288  SUBROUTINE vector_outer(u,v,A)
289  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: u,v
290  REAL(KIND=8), DIMENSION(:,:), INTENT(OUT) :: A
291  INTEGER :: i,j
292 
293  a=0.d0
294  DO i=1,SIZE(u)
295  DO j=1,SIZE(v)
296  a(i,j)=u(i)*v(j)
297  ENDDO
298  ENDDO
299  END SUBROUTINE vector_outer
300 
301 END MODULE util
subroutine, public choldc(a, p)
Definition: util.f90:137
Utility module.
Definition: util.f90:12
subroutine, public init_random_seed()
Random generator initialization routine.
Definition: util.f90:45
integer function, public floordiv(i, j)
Definition: util.f90:241
subroutine, public cdiag(A, d)
Definition: util.f90:230
character(len=40) function, public rstr(x, fm)
Convert a real to string with a given format.
Definition: util.f90:37
integer function lcg(s)
Definition: util.f90:85
subroutine, public printmat(A)
Definition: util.f90:161
subroutine, public diag(A, d)
Definition: util.f90:220
subroutine, public ireduce(A, Ared, n, ind, rind)
Definition: util.f90:275
real(kind=8) function, public mat_trace(A)
Definition: util.f90:111
character(len=20) function, public str(k)
Convert an integer to string.
Definition: util.f90:30
subroutine, public vector_outer(u, v, A)
Definition: util.f90:289
real(kind=8) function, public mat_contract(A, B)
Definition: util.f90:123
subroutine, public reduce(A, Ared, n, ind, rind)
Definition: util.f90:247
subroutine, public triu(A, T)
Definition: util.f90:208
subroutine, public init_one(A)
Initialize a square matrix A as a unit matrix.
Definition: util.f90:100
real(kind=8) function, dimension(size(a, 1), size(a, 2)), public invmat(A)
Definition: util.f90:179
subroutine, public cprintmat(A)
Definition: util.f90:170