-
Notifications
You must be signed in to change notification settings - Fork 184
/
Copy pathstdlib_linalg.fypp
1921 lines (1829 loc) · 84 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
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#:include "common.fypp"
#:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES
#:set RCI_KINDS_TYPES = RC_KINDS_TYPES + INT_KINDS_TYPES
#:set RHS_SUFFIX = ["one","many"]
#:set RHS_SYMBOL = [ranksuffix(r) for r in [1,2]]
#:set RHS_EMPTY = [emptyranksuffix(r) for r in [1,2]]
#:set ALL_RHS = list(zip(RHS_SYMBOL,RHS_SUFFIX,RHS_EMPTY))
#:set EIG_PROBLEM = ["standard", "generalized"]
#:set EIG_FUNCTION = ["geev","ggev"]
#:set EIG_PROBLEM_LIST = list(zip(EIG_PROBLEM, EIG_FUNCTION))
module stdlib_linalg
!!Provides a support for various linear algebra procedures
!! ([Specification](../page/specs/stdlib_linalg.html))
use stdlib_kinds, only: xdp, int8, int16, int32, int64
use stdlib_linalg_constants, only: sp, dp, qp, lk, ilp
use stdlib_error, only: error_stop
use stdlib_optval, only: optval
use stdlib_linalg_state, only: linalg_state_type, linalg_error_handling
implicit none
private
public :: chol
public :: cholesky
public :: det
public :: operator(.det.)
public :: diag
public :: eig
public :: eigh
public :: eigvals
public :: eigvalsh
public :: eye
public :: inv
public :: invert
public :: operator(.inv.)
public :: pinv
public :: pseudoinvert
public :: operator(.pinv.)
public :: lstsq
public :: lstsq_space
public :: norm
public :: mnorm
public :: get_norm
public :: solve
public :: solve_lu
public :: solve_lstsq
public :: trace
public :: svd
public :: svdvals
public :: outer_product
public :: kronecker_product
public :: cross_product
public :: qr
public :: qr_space
public :: schur
public :: schur_space
public :: is_square
public :: is_diagonal
public :: is_symmetric
public :: is_skew_symmetric
public :: hermitian
public :: is_hermitian
public :: is_triangular
public :: is_hessenberg
! Export linalg error handling
public :: linalg_state_type, linalg_error_handling
interface chol
!! version: experimental
!!
!! Computes the Cholesky factorization \( A = L \cdot L^T \), or \( A = U^T \cdot U \).
!! ([Specification](../page/specs/stdlib_linalg.html#chol-compute-the-cholesky-factorization-of-a-rank-2-square-array-matrix))
!!
!!### Summary
!! Pure function interface for computing the Cholesky triangular factors.
!!
!!### Description
!!
!! This interface provides methods for computing the lower- or upper- triangular matrix from the
!! Cholesky factorization of a `real` symmetric or `complex` Hermitian matrix.
!! Supported data types include `real` and `complex`.
!!
!!@note The solution is based on LAPACK's `*POTRF` methods.
!!
#:for rk,rt,ri in RC_KINDS_TYPES
pure module function stdlib_linalg_${ri}$_cholesky_fun(a,lower,other_zeroed) result(c)
!> Input matrix a[m,n]
${rt}$, intent(in) :: a(:,:)
!> [optional] is the lower or upper triangular factor required? Default = lower
logical(lk), optional, intent(in) :: lower
!> [optional] should the unused half of the return matrix be zeroed out? Default: yes
logical(lk), optional, intent(in) :: other_zeroed
!> Output matrix with Cholesky factors c[n,n]
${rt}$ :: c(size(a,1),size(a,2))
end function stdlib_linalg_${ri}$_cholesky_fun
#:endfor
end interface chol
interface cholesky
!! version: experimental
!!
!! Computes the Cholesky factorization \( A = L \cdot L^T \), or \( A = U^T \cdot U \).
!! ([Specification](../page/specs/stdlib_linalg.html#cholesky-compute-the-cholesky-factorization-of-a-rank-2-square-array-matrix))
!!
!!### Summary
!! Pure subroutine interface for computing the Cholesky triangular factors.
!!
!!### Description
!!
!! This interface provides methods for computing the lower- or upper- triangular matrix from the
!! Cholesky factorization of a `real` symmetric or `complex` Hermitian matrix.
!! Supported data types include `real` and `complex`.
!! The factorization is computed in-place if only one matrix argument is present; or returned into
!! a second matrix argument, if present. The `lower` `logical` flag allows to select between upper or
!! lower factorization; the `other_zeroed` optional `logical` flag allows to choose whether the unused
!! part of the triangular matrix should be filled with zeroes.
!!
!!@note The solution is based on LAPACK's `*POTRF` methods.
!!
#:for rk,rt,ri in RC_KINDS_TYPES
pure module subroutine stdlib_linalg_${ri}$_cholesky_inplace(a,lower,other_zeroed,err)
!> Input matrix a[m,n]
${rt}$, intent(inout), target :: a(:,:)
!> [optional] is the lower or upper triangular factor required? Default = lower
logical(lk), optional, intent(in) :: lower
!> [optional] should the unused half of the return matrix be zeroed out? Default: yes
logical(lk), optional, intent(in) :: other_zeroed
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
end subroutine stdlib_linalg_${ri}$_cholesky_inplace
pure module subroutine stdlib_linalg_${ri}$_cholesky(a,c,lower,other_zeroed,err)
!> Input matrix a[n,n]
${rt}$, intent(in) :: a(:,:)
!> Output matrix with Cholesky factors c[n,n]
${rt}$, intent(out) :: c(:,:)
!> [optional] is the lower or upper triangular factor required? Default = lower
logical(lk), optional, intent(in) :: lower
!> [optional] should the unused half of the return matrix be zeroed out? Default: yes
logical(lk), optional, intent(in) :: other_zeroed
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
end subroutine stdlib_linalg_${ri}$_cholesky
#:endfor
end interface cholesky
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
! Identity matrix
interface eye
!! version: experimental
!!
!! Constructs the identity matrix
!! ([Specification](../page/specs/stdlib_linalg.html#eye-construct-the-identity-matrix))
#:for k1, t1 in RCI_KINDS_TYPES
module procedure eye_${t1[0]}$${k1}$
#:endfor
end interface eye
! 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
interface kronecker_product
!! version: experimental
!!
!! Computes the Kronecker product of two arrays of size M1xN1, and of M2xN2, returning an (M1*M2)x(N1*N2) array
!! ([Specification](../page/specs/stdlib_linalg.html#
!! kronecker_product-computes-the-kronecker-product-of-two-matrices))
#:for k1, t1 in RCI_KINDS_TYPES
pure module function kronecker_product_${t1[0]}$${k1}$(A, B) result(C)
${t1}$, intent(in) :: A(:,:), B(:,:)
${t1}$ :: C(size(A,dim=1)*size(B,dim=1),size(A,dim=2)*size(B,dim=2))
end function kronecker_product_${t1[0]}$${k1}$
#:endfor
end interface kronecker_product
! Cross product (of two vectors)
interface cross_product
!! version: experimental
!!
!! Computes the cross product of two vectors, returning a rank-1 and size-3 array
!! ([Specification](../page/specs/stdlib_linalg.html#cross_product-computes-the-cross-product-of-two-3-d-vectors))
#:for k1, t1 in RCI_KINDS_TYPES
pure module function cross_product_${t1[0]}$${k1}$(a, b) result(res)
${t1}$, intent(in) :: a(3), b(3)
${t1}$ :: res(3)
end function cross_product_${t1[0]}$${k1}$
#:endfor
end interface cross_product
! 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
interface hermitian
!! version: experimental
!!
!! Computes the Hermitian version of a rank-2 matrix.
!! For complex matrices, this returns `conjg(transpose(a))`.
!! For real or integer matrices, this returns `transpose(a)`.
!!
!! Usage:
!! ```
!! A = reshape([(1, 2), (3, 4), (5, 6), (7, 8)], [2, 2])
!! AH = hermitian(A)
!! ```
!!
!! [Specification](../page/specs/stdlib_linalg.html#hermitian-compute-the-hermitian-version-of-a-rank-2-matrix)
!!
#:for k1, t1 in RCI_KINDS_TYPES
pure module function hermitian_${t1[0]}$${k1}$(a) result(ah)
${t1}$, intent(in) :: a(:,:)
${t1}$ :: ah(size(a, 2), size(a, 1))
end function hermitian_${t1[0]}$${k1}$
#:endfor
end interface 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
! Solve linear system system Ax=b.
interface solve
!! version: experimental
!!
!! Solves the linear system \( A \cdot x = b \) for the unknown vector \( x \) from a square matrix \( A \).
!! ([Specification](../page/specs/stdlib_linalg.html#solve-solves-a-linear-matrix-equation-or-a-linear-system-of-equations))
!!
!!### Summary
!! Interface for solving a linear system arising from a general matrix.
!!
!!### Description
!!
!! This interface provides methods for computing the solution of a linear matrix system.
!! Supported data types include `real` and `complex`. No assumption is made on the matrix
!! structure.
!! The function can solve simultaneously either one (from a 1-d right-hand-side vector `b(:)`)
!! or several (from a 2-d right-hand-side vector `b(:,:)`) systems.
!!
!!@note The solution is based on LAPACK's generic LU decomposition based solvers `*GESV`.
!!
#:for nd,ndsuf,nde in ALL_RHS
#:for rk,rt,ri in RC_KINDS_TYPES
module function stdlib_linalg_${ri}$_solve_${ndsuf}$(a,b,overwrite_a,err) result(x)
!> Input matrix a[n,n]
${rt}$, intent(inout), target :: a(:,:)
!> Right hand side vector or array, b[n] or b[n,nrhs]
${rt}$, intent(in) :: b${nd}$
!> [optional] Can A data be overwritten and destroyed?
logical(lk), optional, intent(in) :: overwrite_a
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), intent(out) :: err
!> Result array/matrix x[n] or x[n,nrhs]
${rt}$, allocatable, target :: x${nd}$
end function stdlib_linalg_${ri}$_solve_${ndsuf}$
pure module function stdlib_linalg_${ri}$_pure_solve_${ndsuf}$(a,b) result(x)
!> Input matrix a[n,n]
${rt}$, intent(in) :: a(:,:)
!> Right hand side vector or array, b[n] or b[n,nrhs]
${rt}$, intent(in) :: b${nd}$
!> Result array/matrix x[n] or x[n,nrhs]
${rt}$, allocatable, target :: x${nd}$
end function stdlib_linalg_${ri}$_pure_solve_${ndsuf}$
#:endfor
#:endfor
end interface solve
! Solve linear system Ax = b using LU decomposition (subroutine interface).
interface solve_lu
!! version: experimental
!!
!! Solves the linear system \( A \cdot x = b \) for the unknown vector \( x \) from a square matrix \( A \).
!! ([Specification](../page/specs/stdlib_linalg.html#solve-lu-solves-a-linear-matrix-equation-or-a-linear-system-of-equations-subroutine-interface))
!!
!!### Summary
!! Subroutine interface for solving a linear system using LU decomposition.
!!
!!### Description
!!
!! This interface provides methods for computing the solution of a linear matrix system using
!! a subroutine. Supported data types include `real` and `complex`. No assumption is made on the matrix
!! structure. Preallocated space for the solution vector `x` is user-provided, and it may be provided
!! for the array of pivot indices, `pivot`. If all pre-allocated work spaces are provided, no internal
!! memory allocations take place when using this interface.
!! The function can solve simultaneously either one (from a 1-d right-hand-side vector `b(:)`)
!! or several (from a 2-d right-hand-side vector `b(:,:)`) systems.
!!
!!@note The solution is based on LAPACK's generic LU decomposition based solvers `*GESV`.
!!
#:for nd,ndsuf,nde in ALL_RHS
#:for rk,rt,ri in RC_KINDS_TYPES
pure module subroutine stdlib_linalg_${ri}$_solve_lu_${ndsuf}$(a,b,x,pivot,overwrite_a,err)
!> Input matrix a[n,n]
${rt}$, intent(inout), target :: a(:,:)
!> Right hand side vector or array, b[n] or b[n,nrhs]
${rt}$, intent(in) :: b${nd}$
!> Result array/matrix x[n] or x[n,nrhs]
${rt}$, intent(inout), contiguous, target :: x${nd}$
!> [optional] Storage array for the diagonal pivot indices
integer(ilp), optional, intent(inout), target :: pivot(:)
!> [optional] Can A data be overwritten and destroyed?
logical(lk), optional, intent(in) :: overwrite_a
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
end subroutine stdlib_linalg_${ri}$_solve_lu_${ndsuf}$
#:endfor
#:endfor
end interface solve_lu
! Least squares solution to system Ax=b, i.e. such that the 2-norm abs(b-Ax) is minimized.
interface lstsq
!! version: experimental
!!
!! Computes the squares solution to system \( A \cdot x = b \).
!! ([Specification](../page/specs/stdlib_linalg.html#lstsq-computes-the-least-squares-solution-to-a-linear-matrix-equation))
!!
!!### Summary
!! Interface for computing least squares, i.e. the 2-norm \( || (b-A \cdot x ||_2 \) minimizing solution.
!!
!!### Description
!!
!! This interface provides methods for computing the least squares of a linear matrix system.
!! Supported data types include `real` and `complex`.
!!
!!@note The solution is based on LAPACK's singular value decomposition `*GELSD` methods.
!!
#:for nd,ndsuf,nde in ALL_RHS
#:for rk,rt,ri in RC_KINDS_TYPES
module function stdlib_linalg_${ri}$_lstsq_${ndsuf}$(a,b,cond,overwrite_a,rank,err) result(x)
!> Input matrix a[n,n]
${rt}$, intent(inout), target :: a(:,:)
!> Right hand side vector or array, b[n] or b[n,nrhs]
${rt}$, intent(in) :: b${nd}$
!> [optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
real(${rk}$), optional, intent(in) :: cond
!> [optional] Can A,b data be overwritten and destroyed?
logical(lk), optional, intent(in) :: overwrite_a
!> [optional] Return rank of A
integer(ilp), optional, intent(out) :: rank
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
!> Result array/matrix x[n] or x[n,nrhs]
${rt}$, allocatable, target :: x${nd}$
end function stdlib_linalg_${ri}$_lstsq_${ndsuf}$
#:endfor
#:endfor
end interface lstsq
! Least squares solution to system Ax=b, i.e. such that the 2-norm abs(b-Ax) is minimized.
interface solve_lstsq
!! version: experimental
!!
!! Computes the squares solution to system \( A \cdot x = b \).
!! ([Specification](../page/specs/stdlib_linalg.html#solve-lstsq-compute-the-least-squares-solution-to-a-linear-matrix-equation-subroutine-interface))
!!
!!### Summary
!! Subroutine interface for computing least squares, i.e. the 2-norm \( || (b-A \cdot x ||_2 \) minimizing solution.
!!
!!### Description
!!
!! This interface provides methods for computing the least squares of a linear matrix system using
!! a subroutine. Supported data types include `real` and `complex`. If pre-allocated work spaces
!! are provided, no internal memory allocations take place when using this interface.
!!
!!@note The solution is based on LAPACK's singular value decomposition `*GELSD` methods.
!!
#:for nd,ndsuf,nde in ALL_RHS
#:for rk,rt,ri in RC_KINDS_TYPES
module subroutine stdlib_linalg_${ri}$_solve_lstsq_${ndsuf}$(a,b,x,real_storage,int_storage,&
#{if rt.startswith('c')}#cmpl_storage,#{endif}#cond,singvals,overwrite_a,rank,err)
!> Input matrix a[n,n]
${rt}$, intent(inout), target :: a(:,:)
!> Right hand side vector or array, b[n] or b[n,nrhs]
${rt}$, intent(in) :: b${nd}$
!> Result array/matrix x[n] or x[n,nrhs]
${rt}$, intent(inout), contiguous, target :: x${nd}$
!> [optional] real working storage space
real(${rk}$), optional, intent(inout), target :: real_storage(:)
!> [optional] integer working storage space
integer(ilp), optional, intent(inout), target :: int_storage(:)
#:if rt.startswith('complex')
!> [optional] complex working storage space
${rt}$, optional, intent(inout), target :: cmpl_storage(:)
#:endif
!> [optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
real(${rk}$), optional, intent(in) :: cond
!> [optional] list of singular values [min(m,n)], in descending magnitude order, returned by the SVD
real(${rk}$), optional, intent(out), target :: singvals(:)
!> [optional] Can A,b data be overwritten and destroyed?
logical(lk), optional, intent(in) :: overwrite_a
!> [optional] Return rank of A
integer(ilp), optional, intent(out) :: rank
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
end subroutine stdlib_linalg_${ri}$_solve_lstsq_${ndsuf}$
#:endfor
#:endfor
end interface solve_lstsq
! Return the working array space required by the least squares solver
interface lstsq_space
!! version: experimental
!!
!! Computes the integer, real [, complex] working space required by the least-squares solver
!! ([Specification](../page/specs/stdlib_linalg.html#lstsq-space-compute-internal-working-space-requirements-for-the-least-squares-solver))
!!
!!### Description
!!
!! This interface provides sizes of integer, real [, complex] working spaces required by the
!! least-squares solver. These sizes can be used to pre-allocated working arrays in case several
!! repeated least-squares solutions to a same system are sought. If pre-allocated working arrays
!! are provided, no internal allocations will take place.
!!
#:for nd,ndsuf,nde in ALL_RHS
#:for rk,rt,ri in RC_KINDS_TYPES
pure module subroutine stdlib_linalg_${ri}$_lstsq_space_${ndsuf}$(a,b,lrwork,liwork#{if rt.startswith('c')}#,lcwork#{endif}#)
!> Input matrix a[m,n]
${rt}$, intent(in), target :: a(:,:)
!> Right hand side vector or array, b[n] or b[n,nrhs]
${rt}$, intent(in) :: b${nd}$
!> Size of the working space arrays
integer(ilp), intent(out) :: lrwork,liwork#{if rt.startswith('c')}#,lcwork#{endif}#
end subroutine stdlib_linalg_${ri}$_lstsq_space_${ndsuf}$
#:endfor
#:endfor
end interface lstsq_space
! QR factorization of rank-2 array A
interface qr
!! version: experimental
!!
!! Computes the QR factorization of matrix \( A = Q R \).
!! ([Specification](../page/specs/stdlib_linalg.html#qr-compute-the-qr-factorization-of-a-matrix))
!!
!!### Summary
!! Compute the QR factorization of a `real` or `complex` matrix: \( A = Q R \), where \( Q \) is orthonormal
!! and \( R \) is upper-triangular. Matrix \( A \) has size `[m,n]`, with \( m\ge n \).
!!
!!### Description
!!
!! This interface provides methods for computing the QR factorization of a matrix.
!! Supported data types include `real` and `complex`. If a pre-allocated work space
!! is provided, no internal memory allocations take place when using this interface.
!!
!! Given `k = min(m,n)`, one can write \( A = \( Q_1 Q_2 \) \cdot \( \frac{R_1}{0}\) \).
!! The user may want the full problem (provide `shape(Q)==[m,m]`, `shape(R)==[m,n]`) or the reduced
!! problem only: \( A = Q_1 R_1 \) (provide `shape(Q)==[m,k]`, `shape(R)==[k,n]`).
!!
!!@note The solution is based on LAPACK's QR factorization (`*GEQRF`) and ordered matrix output (`*ORGQR`, `*UNGQR`).
!!
#:for rk,rt,ri in RC_KINDS_TYPES
pure module subroutine stdlib_linalg_${ri}$_qr(a,q,r,overwrite_a,storage,err)
!> Input matrix a[m,n]
${rt}$, intent(inout), target :: a(:,:)
!> Orthogonal matrix Q ([m,m], or [m,k] if reduced)
${rt}$, intent(out), contiguous, target :: q(:,:)
!> Upper triangular matrix R ([m,n], or [k,n] if reduced)
${rt}$, intent(out), contiguous, target :: r(:,:)
!> [optional] Can A data be overwritten and destroyed?
logical(lk), optional, intent(in) :: overwrite_a
!> [optional] Provide pre-allocated workspace, size to be checked with qr_space
${rt}$, intent(out), optional, target :: storage(:)
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
end subroutine stdlib_linalg_${ri}$_qr
#:endfor
end interface qr
! Return the working array space required by the QR factorization solver
interface qr_space
!! version: experimental
!!
!! Computes the working array space required by the QR factorization solver
!! ([Specification](../page/specs/stdlib_linalg.html#qr-space-compute-internal-working-space-requirements-for-the-qr-factorization))
!!
!!### Description
!!
!! This interface returns the size of the `real` or `complex` working storage required by the
!! QR factorization solver. The working size only depends on the kind (`real` or `complex`) and size of
!! the matrix being factorized. Storage size can be used to pre-allocate a working array in case several
!! repeated QR factorizations to a same-size matrix are sought. If pre-allocated working arrays
!! are provided, no internal allocations will take place during the factorization.
!!
#:for rk,rt,ri in RC_KINDS_TYPES
pure module subroutine get_qr_${ri}$_workspace(a,lwork,err)
!> Input matrix a[m,n]
${rt}$, intent(in), target :: a(:,:)
!> Minimum workspace size for both operations
integer(ilp), intent(out) :: lwork
!> State return flag. Returns an error if the query failed
type(linalg_state_type), optional, intent(out) :: err
end subroutine get_qr_${ri}$_workspace
#:endfor
end interface qr_space
! Schur decomposition of rank-2 array A
interface schur
!! version: experimental
!!
!! Computes the Schur decomposition of matrix \( A = Z T Z^H \).
!! ([Specification](../page/specs/stdlib_linalg.html#schur-compute-the-schur-decomposition-of-a-matrix))
!!
!!### Summary
!! Compute the Schur decomposition of a `real` or `complex` matrix: \( A = Z T Z^H \), where \( Z \) is
!! orthonormal/unitary and \( T \) is upper-triangular or quasi-upper-triangular. Matrix \( A \) has size `[m,m]`.
!!
!!### Description
!!
!! This interface provides methods for computing the Schur decomposition of a matrix.
!! Supported data types include `real` and `complex`. If a pre-allocated workspace is provided, no internal
!! memory allocations take place when using this interface.
!!
!! The output matrix \( T \) is upper-triangular for `complex` input matrices and quasi-upper-triangular
!! for `real` input matrices (with possible \( 2 \times 2 \) blocks on the diagonal).
!!
!!@note The solution is based on LAPACK's Schur decomposition routines (`*GEES`).
!!
#:for rk,rt,ri in RC_KINDS_TYPES
module subroutine stdlib_linalg_${ri}$_schur(a, t, z, eigvals, overwrite_a, storage, err)
!> Input matrix a[m,m]
${rt}$, intent(inout), target :: a(:,:)
!> Schur form of A: upper-triangular or quasi-upper-triangular matrix T
${rt}$, intent(out), contiguous, target :: t(:,:)
!> Unitary/orthonormal transformation matrix Z
${rt}$, optional, intent(out), contiguous, target :: z(:,:)
!> [optional] Output eigenvalues that appear on the diagonal of T
complex(${rk}$), optional, intent(out), contiguous, target :: eigvals(:)
!> [optional] Can A data be overwritten and destroyed?
logical(lk), optional, intent(in) :: overwrite_a
!> [optional] Provide pre-allocated workspace, size to be checked with schur_space
${rt}$, optional, intent(inout), target :: storage(:)
!> [optional] State return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
end subroutine stdlib_linalg_${ri}$_schur
! Schur decomposition subroutine: real eigenvalue interface
module subroutine stdlib_linalg_real_eig_${ri}$_schur(a,t,z,eigvals,overwrite_a,storage,err)
!> Input matrix a[m,m]
${rt}$, intent(inout), target :: a(:,:)
!> Schur form of A: upper-triangular or quasi-upper-triangular matrix T
${rt}$, intent(out), contiguous, target :: t(:,:)
!> Unitary/orthonormal transformation matrix Z
${rt}$, optional, intent(out), contiguous, target :: z(:,:)
!> Output real eigenvalues that appear on the diagonal of T
real(${rk}$), intent(out), contiguous, target :: eigvals(:)
!> [optional] Provide pre-allocated workspace, size to be checked with schur_space
${rt}$, optional, intent(inout), target :: storage(:)
!> [optional] Can A data be overwritten and destroyed?
logical(lk), optional, intent(in) :: overwrite_a
!> [optional] State return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
end subroutine stdlib_linalg_real_eig_${ri}$_schur
#:endfor
end interface schur
! Return the working array space required by the Schur decomposition solver
interface schur_space
!! version: experimental
!!
!! Computes the working array space required by the Schur decomposition solver
!! ([Specification](../page/specs/stdlib_linalg.html#schur-space-compute-internal-working-space-requirements-for-the-schur-decomposition))
!!
!!### Description
!!
!! This interface returns the size of the `real` or `complex` working storage required by the
!! Schur decomposition solver. The working size only depends on the kind (`real` or `complex`) and size of
!! the matrix being decomposed. Storage size can be used to pre-allocate a working array in case several
!! repeated Schur decompositions of same-size matrices are sought. If pre-allocated working arrays
!! are provided, no internal allocations will take place during the decomposition.
!!
#:for rk,rt,ri in RC_KINDS_TYPES
module subroutine get_schur_${ri}$_workspace(a,lwork,err)
!> Input matrix a[m,m]
${rt}$, intent(in), target :: a(:,:)
!> Minimum workspace size for the decomposition operation
integer(ilp), intent(out) :: lwork
!> State return flag. Returns an error if the query failed
type(linalg_state_type), optional, intent(out) :: err
end subroutine get_schur_${ri}$_workspace
#:endfor
end interface schur_space
interface det
!! version: experimental
!!
!! Computes the determinant of a square matrix
!! ([Specification](../page/specs/stdlib_linalg.html#det-computes-the-determinant-of-a-square-matrix))
!!
!!### Summary
!! Interface for computing matrix determinant.
!!
!!### Description
!!
!! This interface provides methods for computing the determinant of a matrix.
!! Supported data types include `real` and `complex`.
!!
!!@note The provided functions are intended for square matrices only.
!!@note BLAS/LAPACK backends do not currently support extended precision (``xdp``).
!!
!!### Example
!!
!!```fortran
!!
!! real(sp) :: a(3,3), d
!! type(linalg_state_type) :: state
!! a = reshape([1, 2, 3, 4, 5, 6, 7, 8, 9], [3, 3])
!!
!! ! ...
!! d = det(a,err=state)
!! if (state%ok()) then
!! print *, 'Success! det=',d
!! else
!! print *, state%print()
!! endif
!! ! ...
!!```
!!
#:for rk,rt in RC_KINDS_TYPES
#:if rk!="xdp"
module procedure stdlib_linalg_${rt[0]}$${rk}$determinant
module procedure stdlib_linalg_pure_${rt[0]}$${rk}$determinant
#:endif
#:endfor
end interface det
interface operator(.det.)
!! version: experimental
!!
!! Determinant operator of a square matrix
!! ([Specification](../page/specs/stdlib_linalg.html#det-determinant-operator-of-a-square-matrix))
!!
!!### Summary
!! Pure operator interface for computing matrix determinant.
!!
!!### Description
!!
!! This pure operator interface provides a convenient way to compute the determinant of a matrix.
!! Supported data types include real and complex.
!!
!!@note The provided functions are intended for square matrices.
!!@note BLAS/LAPACK backends do not currently support extended precision (``xdp``).
!!
!!### Example
!!
!!```fortran
!!
!! ! ...
!! real(sp) :: matrix(3,3), d
!! matrix = reshape([1, 2, 3, 4, 5, 6, 7, 8, 9], [3, 3])
!! d = .det.matrix
!! ! ...
!!
!!```
!
#:for rk,rt in RC_KINDS_TYPES
#:if rk!="xdp"
module procedure stdlib_linalg_pure_${rt[0]}$${rk}$determinant
#:endif
#:endfor
end interface operator(.det.)
interface
#:for rk,rt in RC_KINDS_TYPES
#:if rk!="xdp"
module function stdlib_linalg_${rt[0]}$${rk}$determinant(a,overwrite_a,err) result(det)
!> Input matrix a[m,n]
${rt}$, intent(inout), target :: a(:,:)
!> [optional] Can A data be overwritten and destroyed?
logical(lk), optional, intent(in) :: overwrite_a
!> State return flag.
type(linalg_state_type), intent(out) :: err
!> Matrix determinant
${rt}$ :: det
end function stdlib_linalg_${rt[0]}$${rk}$determinant
pure module function stdlib_linalg_pure_${rt[0]}$${rk}$determinant(a) result(det)
!> Input matrix a[m,n]
${rt}$, intent(in) :: a(:,:)
!> Matrix determinant
${rt}$ :: det
end function stdlib_linalg_pure_${rt[0]}$${rk}$determinant
#:endif
#:endfor
end interface
! Matrix Inverse: Function interface
interface inv
!! version: experimental
!!
!! Inverse of a square matrix
!! ([Specification](../page/specs/stdlib_linalg.html#inv-inverse-of-a-square-matrix))
!!
!!### Summary
!! This interface provides methods for computing the inverse of a square `real` or `complex` matrix.
!! The inverse \( A^{-1} \) is defined such that \( A \cdot A^{-1} = A^{-1} \cdot A = I_n \).
!!
!!### Description
!!
!! This function interface provides methods that return the inverse of a square matrix.
!! Supported data types include `real` and `complex`.
!! The inverse matrix \( A^{-1} \) is returned as a function result.
!! Exceptions are raised in case of singular matrix or invalid size, and trigger an `error stop` if
!! the state flag `err` is not provided.
!!
!!@note The provided functions are intended for square matrices.
!!
#:for rk,rt,ri in RC_KINDS_TYPES
module function stdlib_linalg_inverse_${ri}$(a,err) result(inva)
!> Input matrix a[n,n]
${rt}$, intent(in) :: a(:,:)
!> Output matrix inverse
${rt}$, allocatable :: inva(:,:)
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
end function stdlib_linalg_inverse_${ri}$
#:endfor
end interface inv
! Matrix Inverse: Subroutine interface - in-place inversion
interface invert
!! version: experimental
!!
!! Inversion of a square matrix
!! ([Specification](../page/specs/stdlib_linalg.html#invert-inversion-of-a-square-matrix))
!!
!!### Summary
!! This interface provides methods for inverting a square `real` or `complex` matrix in-place.
!! The inverse \( A^{-1} \) is defined such that \( A \cdot A^{-1} = A^{-1} \cdot A = I_n \).
!!
!!### Description
!!
!! This subroutine interface provides a way to compute the inverse of a matrix.
!! Supported data types include `real` and `complex`.
!! The user may provide a unique matrix argument `a`. In this case, `a` is replaced by the inverse matrix.
!! on output. Otherwise, one may provide two separate arguments: an input matrix `a` and an output matrix
!! `inva`. In this case, `a` will not be modified, and the inverse is returned in `inva`.
!! Pre-allocated storage may be provided for the array of pivot indices, `pivot`. If all pre-allocated
!! work spaces are provided, no internal memory allocations take place when using this interface.
!!
!!@note The provided subroutines are intended for square matrices.
!!
#:for rk,rt,ri in RC_KINDS_TYPES
module subroutine stdlib_linalg_invert_inplace_${ri}$(a,pivot,err)
!> Input matrix a[n,n]
${rt}$, intent(inout) :: a(:,:)
!> [optional] Storage array for the diagonal pivot indices
integer(ilp), optional, intent(inout), target :: pivot(:)
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
end subroutine stdlib_linalg_invert_inplace_${ri}$
! Compute the square matrix inverse of a
module subroutine stdlib_linalg_invert_split_${ri}$(a,inva,pivot,err)
!> Input matrix a[n,n].
${rt}$, intent(in) :: a(:,:)
!> Inverse matrix a[n,n].
${rt}$, intent(out) :: inva(:,:)
!> [optional] Storage array for the diagonal pivot indices
integer(ilp), optional, intent(inout), target :: pivot(:)
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
end subroutine stdlib_linalg_invert_split_${ri}$
#:endfor
end interface invert
! Matrix Inverse: Operator interface
interface operator(.inv.)
!! version: experimental
!!
!! Inverse operator of a square matrix
!! ([Specification](../page/specs/stdlib_linalg.html#inv-inverse-operator-of-a-square-matrix))
!!
!!### Summary
!! Operator interface for computing the inverse of a square `real` or `complex` matrix.
!!
!!### Description
!!
!! This operator interface provides a convenient way to compute the inverse of a matrix.
!! Supported data types include `real` and `complex`. On input errors or singular matrix,
!! NaNs will be returned.
!!
!!@note The provided functions are intended for square matrices.
!!
#:for rk,rt,ri in RC_KINDS_TYPES
module function stdlib_linalg_inverse_${ri}$_operator(a) result(inva)
!> Input matrix a[n,n]
${rt}$, intent(in) :: a(:,:)
!> Result matrix
${rt}$, allocatable :: inva(:,:)
end function stdlib_linalg_inverse_${ri}$_operator
#:endfor
end interface operator(.inv.)
! Moore-Penrose Pseudo-Inverse: Function interface
interface pinv
!! version: experimental
!!
!! Pseudo-inverse of a matrix
!! ([Specification](../page/specs/stdlib_linalg.html#pinv-moore-penrose-pseudo-inverse-of-a-matrix))
!!
!!### Summary
!! This interface provides methods for computing the Moore-Penrose pseudo-inverse of a matrix.
!! The pseudo-inverse \( A^{+} \) is a generalization of the matrix inverse, computed for square, singular,
!! or rectangular matrices. It is defined such that it satisfies the conditions:
!! - \( A \cdot A^{+} \cdot A = A \)
!! - \( A^{+} \cdot A \cdot A^{+} = A^{+} \)
!! - \( (A \cdot A^{+})^T = A \cdot A^{+} \)
!! - \( (A^{+} \cdot A)^T = A^{+} \cdot A \)
!!
!!### Description
!!
!! This function interface provides methods that return the Moore-Penrose pseudo-inverse of a matrix.
!! Supported data types include `real` and `complex`.
!! The pseudo-inverse \( A^{+} \) is returned as a function result. The computation is based on the
!! singular value decomposition (SVD). An optional relative tolerance `rtol` is provided to control the
!! inclusion of singular values during inversion. Singular values below \( \text{rtol} \cdot \sigma_{\max} \)
!! are treated as zero, where \( \sigma_{\max} \) is the largest singular value. If `rtol` is not provided,
!! a default threshold is applied.
!!
!! Exceptions are raised in case of computational errors or invalid input, and trigger an `error stop`
!! if the state flag `err` is not provided.
!!
!!@note The provided functions are intended for both rectangular and square matrices.
!!
#:for rk,rt,ri in RC_KINDS_TYPES
module function stdlib_linalg_pseudoinverse_${ri}$(a,rtol,err) result(pinva)
!> Input matrix a[m,n]
${rt}$, intent(in), target :: a(:,:)
!> [optional] Relative tolerance for singular value cutoff
real(${rk}$), optional, intent(in) :: rtol
!> [optional] State return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
!> Output matrix pseudo-inverse [n,m]
${rt}$ :: pinva(size(a,2,kind=ilp),size(a,1,kind=ilp))
end function stdlib_linalg_pseudoinverse_${ri}$
#:endfor
end interface pinv
! Moore-Penrose Pseudo-Inverse: Subroutine interface
interface pseudoinvert
!! version: experimental
!!
!! Computation of the Moore-Penrose pseudo-inverse
!! ([Specification](../page/specs/stdlib_linalg.html#pseudoinvert-moore-penrose-pseudo-inverse-of-a-matrix))
!!
!!### Summary
!! This interface provides methods for computing the Moore-Penrose pseudo-inverse of a rectangular
!! or square `real` or `complex` matrix.
!! The pseudo-inverse \( A^{+} \) generalizes the matrix inverse and satisfies the properties:
!! - \( A \cdot A^{+} \cdot A = A \)
!! - \( A^{+} \cdot A \cdot A^{+} = A^{+} \)