-
Notifications
You must be signed in to change notification settings - Fork 184
/
Copy pathstdlib_linalg.fypp
428 lines (382 loc) · 13.7 KB
/
stdlib_linalg.fypp
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
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
#:include "common.fypp"
#:set RCI_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES + INT_KINDS_TYPES
module stdlib_linalg
!!Provides a support for various linear algebra procedures
!! ([Specification](../page/specs/stdlib_linalg.html))
use stdlib_kinds, only: sp, dp, xdp, qp, &
int8, int16, int32, int64
use stdlib_error, only: error_stop
use stdlib_optval, only: optval
implicit none
private
public :: diag
public :: eye
public :: trace
public :: outer_product
public :: repmat
public :: is_square
public :: is_diagonal
public :: is_symmetric
public :: is_skew_symmetric
public :: is_hermitian
public :: is_triangular
public :: is_hessenberg
interface diag
!! version: experimental
!!
!! Creates a diagonal array or extract the diagonal elements of an array
!! ([Specification](../page/specs/stdlib_linalg.html#
!! diag-create-a-diagonal-array-or-extract-the-diagonal-elements-of-an-array))
!
! Vector to matrix
!
#:for k1, t1 in RCI_KINDS_TYPES
module function diag_${t1[0]}$${k1}$(v) result(res)
${t1}$, intent(in) :: v(:)
${t1}$ :: res(size(v),size(v))
end function diag_${t1[0]}$${k1}$
#:endfor
#:for k1, t1 in RCI_KINDS_TYPES
module function diag_${t1[0]}$${k1}$_k(v,k) result(res)
${t1}$, intent(in) :: v(:)
integer, intent(in) :: k
${t1}$ :: res(size(v)+abs(k),size(v)+abs(k))
end function diag_${t1[0]}$${k1}$_k
#:endfor
!
! Matrix to vector
!
#:for k1, t1 in RCI_KINDS_TYPES
module function diag_${t1[0]}$${k1}$_mat(A) result(res)
${t1}$, intent(in) :: A(:,:)
${t1}$ :: res(minval(shape(A)))
end function diag_${t1[0]}$${k1}$_mat
#:endfor
#:for k1, t1 in RCI_KINDS_TYPES
module function diag_${t1[0]}$${k1}$_mat_k(A,k) result(res)
${t1}$, intent(in) :: A(:,:)
integer, intent(in) :: k
${t1}$ :: res(minval(shape(A))-abs(k))
end function diag_${t1[0]}$${k1}$_mat_k
#:endfor
end interface
! Matrix trace
interface trace
!! version: experimental
!!
!! Computes the trace of a matrix
!! ([Specification](../page/specs/stdlib_linalg.html#
!! trace-trace-of-a-matrix))
#:for k1, t1 in RCI_KINDS_TYPES
module procedure trace_${t1[0]}$${k1}$
#:endfor
end interface
! Outer product (of two vectors)
interface outer_product
!! version: experimental
!!
!! Computes the outer product of two vectors, returning a rank-2 array
!! ([Specification](../page/specs/stdlib_linalg.html#
!! outer_product-computes-the-outer-product-of-two-vectors))
#:for k1, t1 in RCI_KINDS_TYPES
pure module function outer_product_${t1[0]}$${k1}$(u, v) result(res)
${t1}$, intent(in) :: u(:), v(:)
${t1}$ :: res(size(u),size(v))
end function outer_product_${t1[0]}$${k1}$
#:endfor
end interface outer_product
! repeat matrix (of 2d array)
interface repmat
!! version: experimental
!!
!! Creates large matrices from a small array, `repmat()` repeats the given values of the array to create the large matrix.
!! ([Specification](../page/specs/stdlib_linalg.html#
!! repmat-creates-large-matrices-from-a-small-array))
#:for k1, t1 in RCI_KINDS_TYPES
pure module function repmat_${t1[0]}$${k1}$(a,m,n) result(res)
${t1}$, intent(in) :: a(:,:)
integer, intent(in) :: m,n
${t1}$ :: res(m*size(a,1),n*size(a,2))
end function repmat_${t1[0]}$${k1}$
#:endfor
end interface repmat
! Check for squareness
interface is_square
!! version: experimental
!!
!! Checks if a matrix (rank-2 array) is square
!! ([Specification](../page/specs/stdlib_linalg.html#
!! is_square-checks-if-a-matrix-is-square))
#:for k1, t1 in RCI_KINDS_TYPES
module procedure is_square_${t1[0]}$${k1}$
#:endfor
end interface is_square
! Check for diagonality
interface is_diagonal
!! version: experimental
!!
!! Checks if a matrix (rank-2 array) is diagonal
!! ([Specification](../page/specs/stdlib_linalg.html#
!! is_diagonal-checks-if-a-matrix-is-diagonal))
#:for k1, t1 in RCI_KINDS_TYPES
module procedure is_diagonal_${t1[0]}$${k1}$
#:endfor
end interface is_diagonal
! Check for symmetry
interface is_symmetric
!! version: experimental
!!
!! Checks if a matrix (rank-2 array) is symmetric
!! ([Specification](../page/specs/stdlib_linalg.html#
!! is_symmetric-checks-if-a-matrix-is-symmetric))
#:for k1, t1 in RCI_KINDS_TYPES
module procedure is_symmetric_${t1[0]}$${k1}$
#:endfor
end interface is_symmetric
! Check for skew-symmetry
interface is_skew_symmetric
!! version: experimental
!!
!! Checks if a matrix (rank-2 array) is skew-symmetric
!! ([Specification](../page/specs/stdlib_linalg.html#
!! is_skew_symmetric-checks-if-a-matrix-is-skew-symmetric))
#:for k1, t1 in RCI_KINDS_TYPES
module procedure is_skew_symmetric_${t1[0]}$${k1}$
#:endfor
end interface is_skew_symmetric
! Check for Hermiticity
interface is_hermitian
!! version: experimental
!!
!! Checks if a matrix (rank-2 array) is Hermitian
!! ([Specification](../page/specs/stdlib_linalg.html#
!! is_hermitian-checks-if-a-matrix-is-hermitian))
#:for k1, t1 in RCI_KINDS_TYPES
module procedure is_hermitian_${t1[0]}$${k1}$
#:endfor
end interface is_hermitian
! Check for triangularity
interface is_triangular
!! version: experimental
!!
!! Checks if a matrix (rank-2 array) is triangular
!! ([Specification](../page/specs/stdlib_linalg.html#
!! is_triangular-checks-if-a-matrix-is-triangular))
#:for k1, t1 in RCI_KINDS_TYPES
module procedure is_triangular_${t1[0]}$${k1}$
#:endfor
end interface is_triangular
! Check for matrix being Hessenberg
interface is_hessenberg
!! version: experimental
!!
!! Checks if a matrix (rank-2 array) is Hessenberg
!! ([Specification](../page/specs/stdlib_linalg.html#
!! is_hessenberg-checks-if-a-matrix-is-hessenberg))
#:for k1, t1 in RCI_KINDS_TYPES
module procedure is_Hessenberg_${t1[0]}$${k1}$
#:endfor
end interface is_hessenberg
contains
!> Version: experimental
!>
!> Constructs the identity matrix.
!> ([Specification](../page/specs/stdlib_linalg.html#eye-construct-the-identity-matrix))
pure function eye(dim1, dim2) result(result)
integer, intent(in) :: dim1
integer, intent(in), optional :: dim2
integer(int8), allocatable :: result(:, :)
integer :: dim2_
integer :: i
dim2_ = optval(dim2, dim1)
allocate(result(dim1, dim2_))
result = 0_int8
do i = 1, min(dim1, dim2_)
result(i, i) = 1_int8
end do
end function eye
#:for k1, t1 in RCI_KINDS_TYPES
function trace_${t1[0]}$${k1}$(A) result(res)
${t1}$, intent(in) :: A(:,:)
${t1}$ :: res
integer :: i
res = 0
do i = 1, minval(shape(A))
res = res + A(i,i)
end do
end function trace_${t1[0]}$${k1}$
#:endfor
#:for k1, t1 in RCI_KINDS_TYPES
pure function is_square_${t1[0]}$${k1}$(A) result(res)
${t1}$, intent(in) :: A(:,:)
logical :: res
res = (size(A,1) == size(A,2))
end function is_square_${t1[0]}$${k1}$
#:endfor
#:for k1, t1 in RCI_KINDS_TYPES
pure function is_diagonal_${t1[0]}$${k1}$(A) result(res)
${t1}$, intent(in) :: A(:,:)
logical :: res
${t1}$, parameter :: zero = 0 !zero of relevant type
integer :: m, n, o, i, j
m = size(A,1)
n = size(A,2)
do j = 1, n !loop over all columns
o = min(j-1,m) !index of row above diagonal (or last row)
do i = 1, o !loop over rows above diagonal
if (A(i,j) /= zero) then
res = .false.
return
end if
end do
do i = o+2, m !loop over rows below diagonal
if (A(i,j) /= zero) then
res = .false.
return
end if
end do
end do
res = .true. !otherwise A is diagonal
end function is_diagonal_${t1[0]}$${k1}$
#:endfor
#:for k1, t1 in RCI_KINDS_TYPES
pure function is_symmetric_${t1[0]}$${k1}$(A) result(res)
${t1}$, intent(in) :: A(:,:)
logical :: res
integer :: n, i, j
if (.not. is_square(A)) then
res = .false.
return !nonsquare matrices cannot be symmetric
end if
n = size(A,1) !symmetric dimension of A
do j = 1, n !loop over all columns
do i = 1, j-1 !loop over all rows above diagonal
if (A(i,j) /= A(j,i)) then
res = .false.
return
end if
end do
end do
res = .true. !otherwise A is symmetric
end function is_symmetric_${t1[0]}$${k1}$
#:endfor
#:for k1, t1 in RCI_KINDS_TYPES
pure function is_skew_symmetric_${t1[0]}$${k1}$(A) result(res)
${t1}$, intent(in) :: A(:,:)
logical :: res
integer :: n, i, j
if (.not. is_square(A)) then
res = .false.
return !nonsquare matrices cannot be skew-symmetric
end if
n = size(A,1) !symmetric dimension of A
do j = 1, n !loop over all columns
do i = 1, j !loop over all rows above diagonal (and diagonal)
if (A(i,j) /= -A(j,i)) then
res = .false.
return
end if
end do
end do
res = .true. !otherwise A is skew-symmetric
end function is_skew_symmetric_${t1[0]}$${k1}$
#:endfor
#:for k1, t1 in (REAL_KINDS_TYPES + INT_KINDS_TYPES)
pure function is_hermitian_${t1[0]}$${k1}$(A) result(res)
${t1}$, intent(in) :: A(:,:)
logical :: res
res = is_symmetric(A) !symmetry and Hermiticity are equivalent for real matrices
end function is_hermitian_${t1[0]}$${k1}$
#:endfor
#:for k1, t1 in CMPLX_KINDS_TYPES
pure function is_hermitian_${t1[0]}$${k1}$(A) result(res)
${t1}$, intent(in) :: A(:,:)
logical :: res
integer :: n, i, j
if (.not. is_square(A)) then
res = .false.
return !nonsquare matrices cannot be Hermitian
end if
n = size(A,1) !symmetric dimension of A
do j = 1, n !loop over all columns
do i = 1, j !loop over all rows above diagonal (and diagonal)
if (A(i,j) /= conjg(A(j,i))) then
res = .false.
return
end if
end do
end do
res = .true. !otherwise A is Hermitian
end function is_hermitian_${t1[0]}$${k1}$
#:endfor
#:for k1, t1 in RCI_KINDS_TYPES
function is_triangular_${t1[0]}$${k1}$(A,uplo) result(res)
${t1}$, intent(in) :: A(:,:)
character, intent(in) :: uplo
logical :: res
${t1}$, parameter :: zero = 0 !zero of relevant type
integer :: m, n, o, i, j
m = size(A,1)
n = size(A,2)
if ((uplo == 'u') .or. (uplo == 'U')) then !check for upper triangularity
do j = 1, n !loop over all columns
o = min(j-1,m) !index of row above diagonal (or last row)
do i = o+2, m !loop over rows below diagonal
if (A(i,j) /= zero) then
res = .false.
return
end if
end do
end do
else if ((uplo == 'l') .or. (uplo == 'L')) then !check for lower triangularity
do j=1,n !loop over all columns
o = min(j-1,m) !index of row above diagonal (or last row)
do i=1,o !loop over rows above diagonal
if (A(i,j) /= zero) then
res = .false.
return
end if
end do
end do
else
call error_stop("ERROR (is_triangular): second argument must be one of {'u','U','l','L'}")
end if
res = .true. !otherwise A is triangular of the requested type
end function is_triangular_${t1[0]}$${k1}$
#:endfor
#:for k1, t1 in RCI_KINDS_TYPES
function is_hessenberg_${t1[0]}$${k1}$(A,uplo) result(res)
${t1}$, intent(in) :: A(:,:)
character, intent(in) :: uplo
logical :: res
${t1}$, parameter :: zero = 0 !zero of relevant type
integer :: m, n, o, i, j
m = size(A,1)
n = size(A,2)
if ((uplo == 'u') .or. (uplo == 'U')) then !check for upper Hessenberg
do j = 1, n !loop over all columns
o = min(j-2,m) !index of row two above diagonal (or last row)
do i = o+4, m !loop over rows two or more below main diagonal
if (A(i,j) /= zero) then
res = .false.
return
end if
end do
end do
else if ((uplo == 'l') .or. (uplo == 'L')) then !check for lower Hessenberg
do j = 1, n !loop over all columns
o = min(j-2,m) !index of row two above diagonal (or last row)
do i = 1, o !loop over rows one or more above main diagonal
if (A(i,j) /= zero) then
res = .false.
return
end if
end do
end do
else
call error_stop("ERROR (is_hessenberg): second argument must be one of {'u','U','l','L'}")
end if
res = .true. !otherwise A is Hessenberg of the requested type
end function is_hessenberg_${t1[0]}$${k1}$
#:endfor
end module stdlib_linalg