The London Perl and Raku Workshop takes place on 26th Oct 2024. If your company depends on Perl, please consider sponsoring and/or attending.

NAME

PDL::GSL::LINALG - PDL interface to linear algebra routines in GSL

SYNOPSIS

use PDL::LiteF;
use PDL::MatrixOps; # for 'x'
use PDL::GSL::LINALG;
my $A = pdl [
  [0.18, 0.60, 0.57, 0.96],
  [0.41, 0.24, 0.99, 0.58],
  [0.14, 0.30, 0.97, 0.66],
  [0.51, 0.13, 0.19, 0.85],
];
my $B = sequence(2,4); # column vectors
LU_decomp(my $lu=$A->copy, my $p=null, my $signum=null);
# transpose so first dim means is vector, higher dims broadcast
LU_solve($lu, $p, $B->transpose, my $x=null);
$x = $x->inplace->transpose; # now can be matrix-multiplied

DESCRIPTION

This is an interface to the linear algebra package present in the GNU Scientific Library. Functions are named as in GSL, but with the initial gsl_linalg_ removed. They are provided in both real and complex double precision.

Currently only LU decomposition interfaces here. Pull requests welcome! #line 60 "LINALG.pm"

FUNCTIONS

LU_decomp

Signature: ([io,phys]A(n,m); indx [o,phys]ipiv(p=CALC($PDL(A)->ndims > 1 ? PDLMIN($PDL(A)->dims[0], $PDL(A)->dims[1]) : 1)); int [o,phys]signum())

LU decomposition of the given (real or complex) matrix.

LU_decomp ignores the bad-value flag of the input ndarrays. It will set the bad-value flag of all output ndarrays if the flag is set for any of the input ndarrays.

LU_solve

Signature: ([phys]LU(n,m); indx [phys]ipiv(p); [phys]B(n); [o,phys]x(n))

Solve A x = B using the LU and permutation from "LU_decomp", real or complex.

LU_solve ignores the bad-value flag of the input ndarrays. It will set the bad-value flag of all output ndarrays if the flag is set for any of the input ndarrays.

LU_det

Signature: ([phys]LU(n,m); int [phys]signum(); [o]det())

Find the determinant from the LU decomp.

LU_det ignores the bad-value flag of the input ndarrays. It will set the bad-value flag of all output ndarrays if the flag is set for any of the input ndarrays.

solve_tridiag

Signature: ([phys]diag(n); [phys]superdiag(n); [phys]subdiag(n); [phys]B(n); [o,phys]x(n))

Solve A x = B where A is a tridiagonal system. Real only, because GSL does not have a complex function.

solve_tridiag ignores the bad-value flag of the input ndarrays. It will set the bad-value flag of all output ndarrays if the flag is set for any of the input ndarrays.

SEE ALSO

PDL

The GSL documentation for linear algebra is online at https://www.gnu.org/software/gsl/doc/html/linalg.html