1LINALG(3) User Contributed Perl Documentation LINALG(3)
2
3
4
6 PDL::GSL::LINALG - PDL interface to linear algebra routines in GSL
7
9 use PDL::LiteF;
10 use PDL::MatrixOps; # for 'x'
11 use PDL::GSL::LINALG;
12 my $A = pdl [
13 [0.18, 0.60, 0.57, 0.96],
14 [0.41, 0.24, 0.99, 0.58],
15 [0.14, 0.30, 0.97, 0.66],
16 [0.51, 0.13, 0.19, 0.85],
17 ];
18 my $B = sequence(2,4); # column vectors
19 LU_decomp(my $lu=$A->copy, my $p=null, my $signum=null);
20 # transpose so first dim means is vector, higher dims broadcast
21 LU_solve($lu, $p, $B->transpose, my $x=null);
22 $x = $x->inplace->transpose; # now can be matrix-multiplied
23
25 This is an interface to the linear algebra package present in the GNU
26 Scientific Library. Functions are named as in GSL, but with the initial
27 "gsl_linalg_" removed. They are provided in both real and complex
28 double precision.
29
30 Currently only LU decomposition interfaces here. Pull requests welcome!
31 #line 59 "LINALG.pm"
32
34 LU_decomp
35 Signature: ([io,phys]A(n,m); indx [o,phys]ipiv(p); int [o,phys]signum())
36
37 LU decomposition of the given (real or complex) matrix.
38
39 LU_decomp ignores the bad-value flag of the input ndarrays. It will
40 set the bad-value flag of all output ndarrays if the flag is set for
41 any of the input ndarrays.
42
43 LU_solve
44 Signature: ([phys]LU(n,m); indx [phys]ipiv(p); [phys]B(n); [o,phys]x(n))
45
46 Solve "A x = B" using the LU and permutation from "LU_decomp", real or
47 complex.
48
49 LU_solve ignores the bad-value flag of the input ndarrays. It will set
50 the bad-value flag of all output ndarrays if the flag is set for any of
51 the input ndarrays.
52
53 LU_det
54 Signature: ([phys]LU(n,m); int [phys]signum(); [o]det())
55
56 Find the determinant from the LU decomp.
57
58 LU_det ignores the bad-value flag of the input ndarrays. It will set
59 the bad-value flag of all output ndarrays if the flag is set for any of
60 the input ndarrays.
61
62 solve_tridiag
63 Signature: ([phys]diag(n); [phys]superdiag(n); [phys]subdiag(n); [phys]B(n); [o,phys]x(n))
64
65 Solve "A x = B" where A is a tridiagonal system. Real only, because GSL
66 does not have a complex function.
67
68 solve_tridiag ignores the bad-value flag of the input ndarrays. It
69 will set the bad-value flag of all output ndarrays if the flag is set
70 for any of the input ndarrays.
71
73 PDL
74
75 The GSL documentation for linear algebra is online at
76 <https://www.gnu.org/software/gsl/doc/html/linalg.html>
77
78
79
80perl v5.36.0 2022-07-22 LINALG(3)