@@ -48,6 +48,7 @@ module stdlib_linalg
48
48
public :: is_diagonal
49
49
public :: is_symmetric
50
50
public :: is_skew_symmetric
51
+ public :: hermitian
51
52
public :: is_hermitian
52
53
public :: is_triangular
53
54
public :: is_hessenberg
@@ -298,6 +299,31 @@ module stdlib_linalg
298
299
#:endfor
299
300
end interface is_hermitian
300
301
302
+ interface hermitian
303
+ !! version: experimental
304
+ !!
305
+ !! Computes the Hermitian version of a rank-2 matrix.
306
+ !! For complex matrices, this returns `conjg(transpose(a))`.
307
+ !! For real or integer matrices, this returns `transpose(a)`.
308
+ !!
309
+ !! Usage:
310
+ !! ```
311
+ !! A = reshape([(1, 2), (3, 4), (5, 6), (7, 8)], [2, 2])
312
+ !! AH = hermitian(A)
313
+ !! ```
314
+ !!
315
+ !! [Specification](../page/specs/stdlib_linalg.html#hermitian)
316
+ !!
317
+
318
+ #:for k1, t1 in RCI_KINDS_TYPES
319
+ module function hermitian_${t1[0]}$${k1}$(a) result(ah)
320
+ ${t1}$, intent(in) :: a(:,:)
321
+ ${t1}$ :: ah(size(a, 2), size(a, 1))
322
+ end function hermitian_${t1[0]}$${k1}$
323
+ #:endfor
324
+
325
+ end interface hermitian
326
+
301
327
302
328
! Check for triangularity
303
329
interface is_triangular
@@ -1471,6 +1497,17 @@ contains
1471
1497
end function is_hermitian_${t1[0]}$${k1}$
1472
1498
#:endfor
1473
1499
1500
+ #:for k1, t1 in RCI_KINDS_TYPES
1501
+ pure module function hermitian_${t1[0]}$${k1}$(a) result(ah)
1502
+ ${t1}$, intent(in) :: a(:,:)
1503
+ ${t1}$ :: ah(size(a, 2), size(a, 1))
1504
+ #:if t1.startswith('complex')
1505
+ ah = conjg(transpose(a))
1506
+ #:else
1507
+ ah = transpose(a)
1508
+ #:endif
1509
+ end function hermitian_${t1[0]}$${k1}$
1510
+ #:endfor
1474
1511
1475
1512
#:for k1, t1 in RCI_KINDS_TYPES
1476
1513
function is_triangular_${t1[0]}$${k1}$(A,uplo) result(res)
@@ -1546,4 +1583,5 @@ contains
1546
1583
end function is_hessenberg_${t1[0]}$${k1}$
1547
1584
#:endfor
1548
1585
1586
+
1549
1587
end module stdlib_linalg
0 commit comments