A Modular Arbitrary-Order Ocean-Atmosphere Model -- Stochastic implementation
Functions/Subroutines
util Module Reference

Utility module. More...

Functions/Subroutines

character(len=20) function, public str (k)
 Convert an integer to string. More...
 
character(len=40) function, public rstr (x, fm)
 Convert a real to string with a given format. More...
 
subroutine, public init_random_seed ()
 Random generator initialization routine. More...
 
subroutine, public init_one (A)
 Initialize a square matrix A as a unit matrix. More...
 
real(kind=8) function, public mat_trace (A)
 
real(kind=8) function, public mat_contract (A, B)
 
subroutine, public choldc (a, p)
 
subroutine, public printmat (A)
 
subroutine, public cprintmat (A)
 
real(kind=8) function, dimension(size(a, 1), size(a, 2)), public invmat (A)
 
subroutine, public triu (A, T)
 
subroutine, public diag (A, d)
 
subroutine, public cdiag (A, d)
 
integer function, public floordiv (i, j)
 
subroutine, public reduce (A, Ared, n, ind, rind)
 
subroutine, public ireduce (A, Ared, n, ind, rind)
 
subroutine, public vector_outer (u, v, A)
 

Detailed Description

Utility module.

Function/Subroutine Documentation

subroutine, public util::cdiag ( complex(kind=16), dimension(:,:), intent(in)  A,
complex(kind=16), dimension(:), intent(out)  d 
)

Definition at line 230 of file util.f90.

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
subroutine, public util::choldc ( real(kind=8), dimension(:,:)  a,
real(kind=8), dimension(:)  p 
)

Definition at line 137 of file util.f90.

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
subroutine, public util::cprintmat ( complex(kind=16), dimension(:,:), intent(in)  A)

Definition at line 170 of file util.f90.

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
subroutine, public util::diag ( real(kind=8), dimension(:,:), intent(in)  A,
real(kind=8), dimension(:), intent(out)  d 
)

Definition at line 220 of file util.f90.

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
integer function, public util::floordiv ( integer  i,
integer  j 
)

Definition at line 241 of file util.f90.

241  INTEGER :: i,j,floordiv
242  floordiv=int(floor(real(i)/real(j)))
243  RETURN
subroutine, public util::init_one ( real(kind=8), dimension(:,:), intent(inout)  A)

Initialize a square matrix A as a unit matrix.

Definition at line 100 of file util.f90.

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 
subroutine, public util::init_random_seed ( )

Random generator initialization routine.

Definition at line 45 of file util.f90.

real(kind=8) function, dimension(size(a,1),size(a,2)), public util::invmat ( real(kind=8), dimension(:,:), intent(in)  A)

Definition at line 179 of file util.f90.

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
subroutine, public util::ireduce ( real(kind=8), dimension(:,:), intent(out)  A,
real(kind=8), dimension(:,:), intent(in)  Ared,
integer, intent(in)  n,
integer, dimension(:), intent(in)  ind,
integer, dimension(:), intent(in)  rind 
)

Definition at line 275 of file util.f90.

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
real(kind=8) function, public util::mat_contract ( real(kind=8), dimension(:,:)  A,
real(kind=8), dimension(:,:)  B 
)

Definition at line 123 of file util.f90.

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
real(kind=8) function, public util::mat_trace ( real(kind=8), dimension(:,:)  A)

Definition at line 111 of file util.f90.

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
subroutine, public util::printmat ( real(kind=8), dimension(:,:), intent(in)  A)

Definition at line 161 of file util.f90.

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
subroutine, public util::reduce ( real(kind=8), dimension(:,:), intent(in)  A,
real(kind=8), dimension(:,:), intent(out)  Ared,
integer, intent(out)  n,
integer, dimension(:), intent(out)  ind,
integer, dimension(:), intent(out)  rind 
)

Definition at line 247 of file util.f90.

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
character(len=40) function, public util::rstr ( real(kind=8), intent(in)  x,
character(len=20), intent(in)  fm 
)

Convert a real to string with a given format.

Definition at line 37 of file util.f90.

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)
character(len=20) function, public util::str ( integer, intent(in)  k)

Convert an integer to string.

Definition at line 30 of file util.f90.

30  INTEGER, INTENT(IN) :: k
31  WRITE (str, *) k
32  str = adjustl(str)
character(len=20) function, public str(k)
Convert an integer to string.
Definition: util.f90:30
subroutine, public util::triu ( real(kind=8), dimension(:,:), intent(in)  A,
real(kind=8), dimension(:,:), intent(out)  T 
)

Definition at line 208 of file util.f90.

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
subroutine, public util::vector_outer ( real(kind=8), dimension(:), intent(in)  u,
real(kind=8), dimension(:), intent(in)  v,
real(kind=8), dimension(:,:), intent(out)  A 
)

Definition at line 289 of file util.f90.

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