File Coverage

lib/PDL/MatrixOps.pd
Criterion Covered Total %
statement 158 203 77.8
branch 72 122 59.0
condition 25 75 33.3
subroutine 9 11 81.8
pod 8 8 100.0
total 272 419 64.9


line stmt bran cond sub pod time code
1             pp_addhdr('
2             #include
3             ');
4              
5             use strict;
6             use warnings;
7             use PDL::Types qw(ppdefs_all);
8              
9             { no warnings 'once'; # pass info back to Makefile.PL
10             $PDL::Core::Dev::EXTRAS{$::PDLMOD}{OBJECT} .= join '', map " $::PDLBASE/$_\$(OBJ_EXT)", qw(eigens simq svd eigen matrix sslib);
11             $PDL::Core::Dev::EXTRAS{$::PDLMOD}{INC} .= qq{ "-I$::PDLBASE"};
12             }
13              
14             pp_addpm({At=>'Top'},<<'EOD');
15             =encoding utf8
16              
17             use strict;
18             use warnings;
19              
20             =head1 NAME
21              
22             PDL::MatrixOps -- Some Useful Matrix Operations
23              
24             =head1 SYNOPSIS
25              
26             $inv = $x->inv;
27              
28             $det = $x->det;
29              
30             ($lu,$perm,$par) = $x->lu_decomp;
31             $y = lu_backsub($lu,$perm,$z); # solve $x x $y = $z
32              
33             =head1 DESCRIPTION
34              
35             PDL::MatrixOps is PDL's built-in matrix manipulation code. It
36             contains utilities for many common matrix operations: inversion,
37             determinant finding, eigenvalue/vector finding, singular value
38             decomposition, etc. PDL::MatrixOps routines are written in a mixture
39             of Perl and C, so that they are reliably present even when there is no
40             external library available (e.g.
41             L or any of the PDL::GSL family of modules).
42              
43             Matrix manipulation, particularly with large matrices, is a
44             challenging field and no one algorithm is suitable in all cases. The
45             utilities here use general-purpose algorithms that work acceptably for
46             many cases but might not scale well to very large or pathological
47             (near-singular) matrices.
48              
49             Except as noted, the matrices are PDLs whose 0th dimension ranges over
50             column and whose 1st dimension ranges over row. The matrices appear
51             correctly when printed.
52              
53             These routines should work OK with L objects
54             as well as with normal PDLs.
55              
56             =head1 TIPS ON MATRIX OPERATIONS
57              
58             Like most computer languages, PDL addresses matrices in (column,row)
59             order in most cases; this corresponds to (X,Y) coordinates in the
60             matrix itself, counting rightwards and downwards from the upper left
61             corner. This means that if you print a PDL that contains a matrix,
62             the matrix appears correctly on the screen, but if you index a matrix
63             element, you use the indices in the reverse order that you would in a
64             math textbook. If you prefer your matrices indexed in (row, column)
65             order, you can try using the L object, which
66             includes an implicit exchange of the first two dimensions but should
67             be compatible with most of these matrix operations. TIMTOWDTI.)
68              
69             Matrices, row vectors, and column vectors can be multiplied with the 'x'
70             operator (which is, of course, broadcastable):
71              
72             $m3 = $m1 x $m2;
73             $col_vec2 = $m1 x $col_vec1;
74             $row_vec2 = $row_vec1 x $m1;
75             $scalar = $row_vec x $col_vec;
76              
77             Because of the (column,row) addressing order, 1-D PDLs are treated as
78             _row_ vectors; if you want a _column_ vector you must add a dummy dimension:
79              
80             $rowvec = pdl(1,2); # row vector
81             $colvec = $rowvec->slice('*1'); # 1x2 column vector
82             $matrix = pdl([[3,4],[6,2]]); # 2x2 matrix
83             $rowvec2 = $rowvec x $matrix; # right-multiplication by matrix
84             $colvec = $matrix x $colvec; # left-multiplication by matrix
85             $m2 = $matrix x $rowvec; # Throws an error
86              
87             Implicit broadcasting works correctly with most matrix operations, but
88             you must be extra careful that you understand the dimensionality. In
89             particular, matrix multiplication and other matrix ops need nx1 PDLs
90             as row vectors and 1xn PDLs as column vectors. In most cases you must
91             explicitly include the trailing 'x1' dimension in order to get the expected
92             results when you broadcast over multiple row vectors.
93              
94             When broadcasting over matrices, it's very easy to get confused about
95             which dimension goes where. It is useful to include comments with
96             every expression, explaining what you think each dimension means:
97              
98             $x = xvals(360)*3.14159/180; # (angle)
99             $rot = cat(cat(cos($x),sin($x)), # rotmat: (col,row,angle)
100             cat(-sin($x),cos($x)));
101              
102             =head1 ACKNOWLEDGEMENTS
103              
104             MatrixOps includes algorithms and pre-existing code from several
105             origins. In particular, C is the work of Stephen Moshier,
106             C uses an SVD subroutine written by Bryant Marks, and C
107             uses a subset of the Small Scientific Library by Kenneth Geisshirt.
108             They are free software, distributable under same terms as PDL itself.
109              
110              
111             =head1 NOTES
112              
113             This is intended as a general-purpose linear algebra package for
114             small-to-mid sized matrices. The algorithms may not scale well to
115             large matrices (hundreds by hundreds) or to near singular matrices.
116              
117             If there is something you want that is not here, please add and
118             document it!
119              
120             =cut
121 70     70   621  
  70         158  
  70         5349  
122 70     70   397 use Carp;
  70     0   146  
  70         287895  
123             use strict;
124              
125             EOD
126              
127              
128             ######################################################################
129              
130             pp_add_exported('','identity');
131             pp_addpm(<<'EOD');
132             =head2 identity
133              
134             =for sig
135              
136             Signature: (n; [o]a(n,n))
137              
138             =for ref
139              
140             Return an identity matrix of the specified size. If you hand in a
141             scalar, its value is the size of the identity matrix; if you hand in a
142             dimensioned PDL, the 0th dimension is the first two dimensions of the
143             matrix, with higher dimensions preserved.
144              
145             =cut
146              
147             sub identity {
148 0     65 1 0 my $n = shift;
149 65 100       4703 my $out =
    100          
150             !(my $was_pdl = UNIVERSAL::isa($n,'PDL')) ? zeroes($n,$n) :
151             $n->getndims == 0 ? zeroes($n->type, $n->at(0),$n->at(0)) :
152             undef;
153 65 100       659 if (!defined $out) {
154 65         188 my @dims = $n->dims;
155 50         186 $out = zeroes($n->type, @dims[0, 0, 2..$#dims]);
156             }
157 50         151 $out->diagonal(0,1)++;
158 65 100       237 $was_pdl ? bless $out, ref($n) : $out;
159             }
160             EOD
161              
162             ######################################################################
163             pp_add_exported('','stretcher');
164             pp_addpm(<<'EOD');
165              
166             =head2 stretcher
167              
168             =for sig
169              
170             Signature: (a(n); [o]b(n,n))
171              
172             =for usage
173              
174             $mat = stretcher($eigenvalues);
175              
176             =for ref
177              
178             Return a diagonal matrix with the specified diagonal elements.
179             Preserves higher dimensions.
180             As of 2.096, it will also have the same datatype.
181              
182             =cut
183              
184 65     4 1 581 sub stretcher {
185 4         8 my $in = shift;
186 4         12 my $out = zeroes($in->type, $in->dim(0), $in->dims);
187 4         17 $out->diagonal(0,1) += $in;
188             $out;
189             }
190              
191             EOD
192              
193             ######################################################################
194              
195             pp_add_exported('','inv');
196             pp_addpm(<<'EOD');
197              
198             =head2 inv
199              
200             =for sig
201              
202             Signature: (a(m,m); sv opt )
203              
204             =for usage
205              
206             $a1 = inv($a, {$opt});
207              
208             =for ref
209              
210             Invert a square matrix.
211              
212             You feed in an NxN matrix in $a, and get back its inverse (if it
213             exists). The code is inplace-aware, so you can get back the inverse
214             in $a itself if you want -- though temporary storage is used either
215             way. You can cache the LU decomposition in an output option variable.
216              
217             C uses C by default; that is a numerically stable
218             (pivoting) LU decomposition method.
219              
220              
221             OPTIONS:
222              
223             =over 3
224              
225             =item * s
226              
227             Boolean value indicating whether to complain if the matrix is singular. If
228             this is false, singular matrices cause inverse to barf. If it is true, then
229             singular matrices cause inverse to return undef.
230              
231             =item * lu (I/O)
232              
233             This value contains a list ref with the LU decomposition, permutation,
234             and parity values for C<$a>. If you do not mention the key, or if the
235             value is undef, then inverse calls C. If the key exists with
236             an undef value, then the output of C is stashed here (unless
237             the matrix is singular). If the value exists, then it is assumed to
238             hold the LU decomposition.
239              
240             =item * det (Output)
241              
242             If this key exists, then the determinant of C<$a> get stored here,
243             whether or not the matrix is singular.
244              
245             =back
246              
247             =cut
248              
249 4     47 1 29 *PDL::inv = \&inv;
250 47         576 sub inv {
251 47 100       104 my $x = shift;
252             my $opt = shift;
253 47 50 33     110 $opt = {} unless defined($opt);
      33        
254              
255             barf "inverse needs a square PDL as a matrix\n"
256             unless(UNIVERSAL::isa($x,'PDL') &&
257             $x->dims >= 2 &&
258             $x->dim(0) == $x->dim(1)
259 47         571 );
260 47 100 66     92  
      100        
261             my ($lu,$perm,$par);
262             if(exists($opt->{lu}) &&
263 47         353 ref $opt->{lu} eq 'ARRAY' &&
  37         47  
264             ref $opt->{lu}->[0] eq 'PDL') {
265 37         115 ($lu,$perm,$par) = @{$opt->{lu}};
266 10         67 } else {
267 10 100       39 ($lu,$perm,$par) = lu_decomp($x);
268             @{$opt->{lu}} = ($lu,$perm,$par)
269             if(ref $opt->{lu} eq 'ARRAY');
270 1 50       5 }
271              
272 47 50       205 my $det = (defined $lu) ? $lu->diagonal(0,1)->prodover * $par : pdl(0);
273             $opt->{det} = $det
274 47 100 100     651 if exists($opt->{det});
275              
276 47 50       350 unless($det->nelem > 1 || $det) {
277 1         16 return undef
278             if $opt->{s};
279             barf("PDL::inv: got a singular matrix or LU decomposition\n");
280 0         0 }
281              
282 46 100       608 my $out = lu_backsub($lu,$perm,$par,identity($x))->transpose->sever;
283              
284             return $out
285 46         569 unless($x->is_inplace);
286 3         14  
287             $x .= $out;
288             $x;
289             }
290              
291             EOD
292              
293             ######################################################################
294              
295             pp_add_exported('','det');
296             pp_addpm(<<'EOD');
297              
298             =head2 det
299              
300             =for sig
301              
302             Signature: (a(m,m); sv opt)
303              
304             =for usage
305              
306             $det = det($a,{opt});
307              
308             =for ref
309              
310             Determinant of a square matrix using LU decomposition (for large matrices)
311              
312             You feed in a square matrix, you get back the determinant. Some
313             options exist that allow you to cache the LU decomposition of the
314             matrix (note that the LU decomposition is invalid if the determinant
315             is zero!). The LU decomposition is cacheable, in case you want to
316             re-use it. This method of determinant finding is more rapid than
317             recursive-descent on large matrices, and if you reuse the LU
318             decomposition it's essentially free.
319              
320             OPTIONS:
321              
322             =over 3
323              
324             =item * lu (I/O)
325              
326             Provides a cache for the LU decomposition of the matrix. If you
327             provide the key but leave the value undefined, then the LU decomposition
328             goes in here; if you put an LU decomposition here, it will be used and
329             the matrix will not be decomposed again.
330              
331             =back
332              
333             =cut
334              
335             *PDL::det = \&det;
336 3     42 1 28 sub det {
337 42 100       133 my ($x, $opt) = @_;
338 42         102 $opt = {} unless defined($opt);
339 42 100 100     77 my($lu,$perm,$par);
340 42         346 if(exists ($opt->{lu}) and (ref $opt->{lu} eq 'ARRAY')) {
  1         2  
341             ($lu,$perm,$par) = @{$opt->{lu}};
342 1         2 } else {
343             ($lu,$perm,$par) = lu_decomp($x);
344 41 100       141 $opt->{lu} = [$lu,$perm,$par]
345             if(exists($opt->{lu}));
346 41 100       214 }
347             defined $lu ? $lu->diagonal(0,1)->prodover * $par : PDL->zeroes(sbyte,1);
348             }
349              
350             EOD
351              
352             ######################################################################
353              
354             pp_add_exported('','determinant');
355             pp_addpm(<<'EOD');
356              
357             =head2 determinant
358              
359             =for sig
360              
361             Signature: (a(m,m))
362              
363             =for usage
364              
365             $det = determinant($x);
366              
367             =for ref
368              
369             Determinant of a square matrix, using recursive descent (broadcastable).
370              
371             This is the traditional, robust recursive determinant method taught in
372             most linear algebra courses. It scales like C (and hence is
373             pitifully slow for large matrices) but is very robust because no
374             division is involved (hence no division-by-zero errors for singular
375             matrices). It's also broadcastable, so you can find the determinants of
376             a large collection of matrices all at once if you want.
377              
378             Matrices up to 3x3 are handled by direct multiplication; larger matrices
379             are handled by recursive descent to the 3x3 case.
380              
381             The LU-decomposition method L is faster in isolation for
382             single matrices larger than about 4x4, and is much faster if you end up
383             reusing the LU decomposition of C<$a> (NOTE: check performance and
384             broadcasting benchmarks with new code).
385              
386             =cut
387              
388             *PDL::determinant = \&determinant;
389 42     6 1 184 sub determinant {
390 6         18 my($x) = shift;
391             my($n);
392 6 50 33     7 return undef unless(
      33        
393             UNIVERSAL::isa($x,'PDL') &&
394             $x->getndims >= 2 &&
395             ($n = $x->dim(0)) == $x->dim(1)
396             );
397 6 50       68  
398 6 50       22 return $x->clump(2) if($n==1);
399 6         13 if($n==2) {
400 0         0 my($y) = $x->clump(2);
401             return $y->index(0)*$y->index(3) - $y->index(1)*$y->index(2);
402 0 100       0 }
403 6         12 if($n==3) {
404 5         15 my($y) = $x->clump(2);
405 5         47 my $y3 = $y->index(3);
406 5         41 my $y4 = $y->index(4);
407 5         31 my $y5 = $y->index(5);
408 5         32 my $y6 = $y->index(6);
409 5         31 my $y7 = $y->index(7);
410             my $y8 = $y->index(8);
411              
412 5         75 return (
413             $y->index(0) * ( $y4 * $y8 - $y5 * $y7 )
414             + $y->index(1) * ( $y5 * $y6 - $y3 * $y8 )
415             + $y->index(2) * ( $y3 * $y7 - $y4 * $y6 )
416             );
417             }
418 5         106  
419 1         1 my($i);
420             my($sum) = zeroes($x->slice('(0),(0)'));
421              
422 1         4 # Do middle submatrices
423 1         7 for $i(1..$n-2) {
424 2 50       9 my $el = $x->slice("($i),(0)");
425 2         6 next if( ($el==0)->all ); # Optimize away unnecessary recursion
426             $sum += $el * (1-2*($i%2)) *
427             determinant($x->slice("0:".($i-1).",1:-1")->
428             append($x->slice(($i+1).":-1,1:-1")));
429             }
430              
431 2         13 # Do beginning and end submatrices
432 1         5 $sum += $x->slice("(0),(0)") * determinant($x->slice('1:-1,1:-1'));
433             $sum -= $x->slice("(-1),(0)") * determinant($x->slice('0:-2,1:-1')) * (1 - 2*($n % 2));
434 1         23  
435             return $sum;
436             }
437              
438             EOD
439              
440             my $TRIANGULAR_REDODIMS = <<'EOF';
441             float n = (sqrtf(1 + 8*$SIZE(m)) - 1)/2;
442             $SIZE(n) = roundf(n);
443             if (fabsf($SIZE(n) - n) > 0.0001)
444             $CROAK("Non-triangular vector size=%"IND_FLAG, $SIZE(m));
445             EOF
446              
447             pp_def("eigens_sym",
448             HandleBad => 0,
449             Pars => '[phys]a(m); [o,phys]ev(n,n); [o,phys]e(n)',
450             RedoDimsCode => $TRIANGULAR_REDODIMS,
451             GenericTypes => ['D'],
452             Code => <<'EOF',
453             extern void eigens( double *A, double *EV, double *E, int N );
454             eigens($P(a), $P(ev), $P(e), $SIZE(n));
455             EOF
456             PMCode => <<'EOF',
457             sub PDL::eigens_sym {
458             my ($x) = @_;
459             my @d = $x->dims;
460             barf "Need real square matrix for eigens_sym"
461             if @d < 2 or $d[0] != $d[1];
462             my ($sym) = 0.5*($x + $x->transpose);
463             my ($err) = PDL::max(abs($sym));
464             barf "Need symmetric component non-zero for eigens_sym" if $err == 0;
465             $err = PDL::max(abs($x-$sym))/$err;
466             warn "Using symmetrized version of the matrix in eigens_sym"
467             if $err > 1e-5 && $PDL::debug;
468             PDL::_eigens_sym_int($sym->squaretotri, my $ev=PDL->null, my $e=PDL->null);
469             return ($ev->transpose, $e) if wantarray;
470             $e; #just eigenvalues
471             }
472             EOF
473             Doc => <<'EOF',
474             =for ref
475              
476             Eigenvalues and -vectors of a symmetric square matrix. If passed
477             an asymmetric matrix, the routine will warn and symmetrize it, by taking
478             the average value. That is, it will solve for 0.5*($a+$a->transpose).
479              
480             It's broadcastable, so if C<$a> is 3x3x100, it's treated as 100 separate 3x3
481             matrices, and both C<$ev> and C<$e> get extra dimensions accordingly.
482              
483             If called in scalar context it hands back only the eigenvalues. Ultimately,
484             it should switch to a faster algorithm in this case (as discarding the
485             eigenvectors is wasteful).
486              
487             The algorithm used is due to J. von Neumann, which was a rediscovery of
488             L .
489              
490             The eigenvectors are returned in COLUMNS of the returned PDL. That
491             makes it slightly easier to access individual eigenvectors, since the
492             0th dim of the output PDL runs across the eigenvectors and the 1st dim
493             runs across their components.
494              
495             ($ev,$e) = eigens_sym $x; # Make eigenvector matrix
496             $vector = $ev->slice($n); # Select nth eigenvector as a column-vector
497             $vector = $ev->slice("($n)"); # Select nth eigenvector as a row-vector
498              
499             As of 2.096, the eigenvalues are returned in ascending order.
500              
501             To compare with L:
502              
503             use PDL::LinearAlgebra;
504             ($val, $rvec) = msymeigen($A = pdl([3,4], [4,-3]),1,1);
505             print $val->slice(1) * $rvec->slice(1);
506             #[
507             # [ -4.472136]
508             # [ -2.236068]
509             #]
510             print $A x $rvec->slice(1);
511             #[
512             # [ -4.472136]
513             # [ -2.236068]
514             #]
515             ($rvec, $val) = eigens_sym($A); # note return values other way round
516             # otherwise the same
517              
518             =for usage
519              
520             ($ev, $e) = eigens_sym($x); # e-vects & e-values
521             $e = eigens_sym($x); # just eigenvalues
522             EOF
523             );
524              
525             ######################################################################
526             ### eigens
527             ###
528             pp_def("eigens",
529             HandleBad => 0,
530             Pars => '[phys]a(n,n); complex [o,phys]ev(n,n); complex [o,phys]e(n)',
531             GenericTypes => ['D'],
532             CHeader => '#include "eigen.h"',
533             Code => <<'EOF',
534             char *ret = Eigen($SIZE(n), $P(a), 20*$SIZE(n), 1e-13, $P(e), $P(ev));
535             if (ret) $CROAK("%s", ret);
536             EOF
537             PMCode => <<'EOF',
538             sub PDL::eigens {
539             my ($x) = @_;
540             my @d = $x->dims;
541             barf "Need real square matrix for eigens"
542             if @d < 2 or $d[0] != $d[1];
543             my $deviation = PDL::max(abs($x - $x->transpose))/PDL::max(abs($x));
544             if ( $deviation <= 1e-5 ) {
545             PDL::_eigens_sym_int($x->squaretotri, my $ev=PDL->null, my $e=PDL->null);
546             return $ev->transpose, $e if wantarray;
547             return $e; #just eigenvalues
548             } else {
549             PDL::_eigens_int($x, my $ev=PDL->null, my $e=PDL->null);
550             return $ev, $e if wantarray;
551             return $e; #just eigenvalues
552             }
553             }
554             EOF
555             Doc => <<'EOF',
556             =for ref
557              
558             Complex eigenvalues and -vectors of a real square matrix.
559              
560             (See also L<"eigens_sym"|/eigens_sym>, for eigenvalues and -vectors
561             of a real, symmetric, square matrix).
562              
563             The eigens function will attempt to compute the eigenvalues and
564             eigenvectors of a square matrix with real components. If the matrix
565             is symmetric, the same underlying code as L<"eigens_sym"|/eigens_sym>
566             is used. If asymmetric, the eigenvalues and eigenvectors are computed
567             with algorithms from the sslib library. These are a slightly modified
568             version of EISPACK's C code.
569              
570             Not all square matrices are diagonalizable. If you feed in a
571             non-diagonalizable matrix, then the algorithm may fail, in which case
572             an exception will be thrown.
573              
574             C is broadcastable, so you can solve 100 eigenproblems by
575             feeding in a 3x3x100 array. Both C<$ev> and C<$e> get extra dimensions accordingly.
576              
577             If called in scalar context C hands back only the eigenvalues. This
578             is somewhat wasteful, as it calculates the eigenvectors anyway.
579              
580             The eigenvectors are returned in COLUMNS of the returned PDL (ie the
581             the 0 dimension). That makes it slightly easier to access individual
582             eigenvectors, since the 0th dim of the output PDL runs across the
583             eigenvectors and the 1st dim runs across their components.
584              
585             ($ev,$e) = eigens $x; # Make eigenvector matrix
586             $vector = $ev->slice($n); # Select nth eigenvector as a column-vector
587             $vector = $ev->slice("($n)"); # Select nth eigenvector as a row-vector
588              
589             To compare with L:
590              
591             use PDL::LinearAlgebra;
592             ($val, $lvec, $rvec) = meigen($A = pdl([4,-1], [2,1]),1,1);
593             print $val->slice(1) * $rvec->slice(1);
594             #[
595             # [0.894427190999916]
596             # [ 1.78885438199983]
597             #]
598             print $A x $rvec->slice(1);
599             #[
600             # [0.894427190999916]
601             # [ 1.78885438199983]
602             #]
603             ($rvec, $val) = eigens($A); # note return values other way round
604             # otherwise the same
605              
606             =for usage
607              
608             ($ev, $e) = eigens($x); # e'vects & e'vals
609             $e = eigens($x); # just eigenvalues
610             EOF
611             );
612              
613             ######################################################################
614             ### svd
615             pp_def(
616             "svd",
617             HandleBad => 0,
618             Pars => 'a(n,m); [t]w(wsize=CALC($SIZE(n) * ($SIZE(m) + $SIZE(n)))); [o]u(n,m); [o,phys]z(n); [o]v(n,n);',
619             GenericTypes => ['D'],
620             RedoDimsCode => '
621             if ($SIZE(m)<$SIZE(n))
622             $CROAK("svd requires input ndarrays to have m >= n; you have m=%td and n=%td. Try inputting the transpose. See the docs for svd.",$SIZE(m),$SIZE(n));
623             ',
624             Code => '
625             extern void SVD( double *W, double *Z, int nRow, int nCol );
626             double *t = $P(w);
627             loop (m,n) %{
628             *t++ = $a();
629             %}
630             SVD($P(w), $P(z), $SIZE(m), $SIZE(n));
631             loop (n) %{
632             $z() = sqrt($z());
633             %}
634             t = $P(w);
635             loop (m,n) %{
636             $u() = *t++/$z();
637             %}
638             loop (n1,n0) %{
639             $v() = *t++;
640             %}
641             ',
642             , Doc => q{
643             =for ref
644              
645             Singular value decomposition of a matrix.
646              
647             C is broadcastable.
648              
649             Given an m x n matrix C<$a> that has m rows and n columns (m >= n),
650             C computes matrices C<$u> and C<$v>, and a vector of the singular
651             values C<$s>. Like most implementations, C computes what is
652             commonly referred to as the "thin SVD" of C<$a>, such that C<$u> is m
653             x n, C<$v> is n x n, and there are <=n singular values. As long as m
654             >= n, the original matrix can be reconstructed as follows:
655              
656             ($u,$s,$v) = svd($x);
657             $ess = zeroes($x->dim(0),$x->dim(0));
658             $ess->slice("$_","$_").=$s->slice("$_") foreach (0..$x->dim(0)-1); #generic diagonal
659             $a_copy = $u x $ess x $v->transpose;
660              
661             If m==n, C<$u> and C<$v> can be thought of as rotation matrices that
662             convert from the original matrix's singular coordinates to final
663             coordinates, and from original coordinates to singular coordinates,
664             respectively, and $ess is a diagonal scaling matrix.
665              
666             If n>m, C will barf. This can be avoided by passing in the
667             transpose of C<$a>, and reconstructing the original matrix like so:
668              
669             ($u,$s,$v) = svd($x->transpose);
670             $ess = zeroes($x->dim(1),$x->dim(1));
671             $ess->slice($_,$_).=$s->slice($_) foreach (0..$x->dim(1)-1); #generic diagonal
672             $x_copy = $v x $ess x $u->transpose;
673              
674             EXAMPLE
675              
676             The computing literature has loads of examples of how to use SVD.
677             Here's a trivial example (used in L)
678             of how to make a matrix less, er, singular, without changing the
679             orientation of the ellipsoid of transformation:
680              
681             { my($r1,$s,$r2) = svd $x;
682             $s++; # fatten all singular values
683             $r2 *= $s; # implicit broadcasting for cheap mult.
684             $x .= $r2 x $r1; # a gets r2 x ess x r1
685             }
686             },);
687              
688             ######################################################################
689             pp_add_exported('','lu_decomp');
690             pp_addpm(<<'EOD');
691              
692             =head2 lu_decomp
693              
694             =for sig
695              
696             Signature: (a(m,m); [o]lu(m,m); [o]perm(m); [o]parity)
697              
698             =for ref
699              
700             LU decompose a matrix, with row permutation
701              
702             =for usage
703              
704             ($lu, $perm, $parity) = lu_decomp($x);
705              
706             $lu = lu_decomp($x, $perm, $par); # $perm and $par are outputs!
707              
708             lu_decomp($x->inplace,$perm,$par); # Everything in place.
709              
710             =for description
711              
712             C returns an LU decomposition of a square matrix,
713             using Crout's method with partial pivoting. It's ported
714             from I. The partial pivoting keeps it
715             numerically stable but means a little more overhead from
716             broadcasting.
717              
718             C decomposes the input matrix into matrices L and
719             U such that LU = A, L is a subdiagonal matrix, and U is a
720             superdiagonal matrix. By convention, the diagonal of L is
721             all 1's.
722              
723             The single output matrix contains all the variable elements
724             of both the L and U matrices, stacked together. Because the
725             method uses pivoting (rearranging the lower part of the
726             matrix for better numerical stability), you have to permute
727             input vectors before applying the L and U matrices. The
728             permutation is returned either in the second argument or, in
729             list context, as the second element of the list. You need
730             the permutation for the output to make any sense, so be sure
731             to get it one way or the other.
732              
733             LU decomposition is the answer to a lot of matrix questions,
734             including inversion and determinant-finding, and C
735             is used by L.
736              
737             If you pass in C<$perm> and C<$parity>, they either must be
738             predeclared PDLs of the correct size ($perm is an n-vector,
739             C<$parity> is a scalar) or scalars.
740              
741             If the matrix is singular, then the LU decomposition might
742             not be defined; in those cases, C silently returns
743             undef. Some singular matrices LU-decompose just fine, and
744             those are handled OK but give a zero determinant (and hence
745             can't be inverted).
746              
747             C uses pivoting, which rearranges the values in the
748             matrix for more numerical stability. This makes it really
749             good for large and even near-singular matrices. There is
750             a non-pivoting version C available which is
751             from 5 to 60 percent faster for typical problems at
752             the expense of failing to compute a result in some cases.
753              
754             Now that the C is broadcasted, it is the recommended
755             LU decomposition routine. It no longer falls back to C.
756              
757             C is ported from I to PDL. It
758             should probably be implemented in C.
759              
760             =cut
761              
762             *PDL::lu_decomp = \&lu_decomp;
763              
764 1     58 1 12 sub lu_decomp {
765 58         480 my($in) = shift;
766 58         101 my($permute) = shift;
767 58         92 my($parity) = shift;
768             my($sing_ok) = shift;
769 58         92  
770             my $TINY = 1e-30;
771 58 50 33     94  
      33        
772             barf("lu_decomp requires a square (2D) PDL\n")
773             if(!UNIVERSAL::isa($in,'PDL') ||
774             $in->ndims < 2 ||
775             $in->dim(0) != $in->dim(1));
776 58         788  
777 58         164 my($n) = $in->dim(0);
  58         112  
778             my($n1) = $n; $n1--;
779 58         84  
780 58 100       207 my($inplace) = $in->is_inplace;
781             my($out) = ($inplace) ? $in : $in->copy;
782 58 50       257  
783 58 0 0     170 if(defined $permute) {
      0        
784             barf('lu_decomp: permutation vector must match the matrix')
785             if(!UNIVERSAL::isa($permute,'PDL') ||
786             $permute->ndims != 1 ||
787 0         0 $permute->dim(0) != $out->dim(0));
788             $permute .= PDL->xvals($in->dim(0));
789 0         0 } else {
790             $permute = $in->slice("(0)")->xvals;
791             }
792 58 50       244  
793 58 0 0     437 if(defined $parity) {
794             barf('lu_decomp: parity must be a scalar PDL')
795             if(!UNIVERSAL::isa($parity,'PDL') ||
796 0         0 $parity->dim(0) != 1);
797             $parity .= 1.0;
798 0         0 } else {
799             $parity = $in->slice('(0),(0)')->ones;
800             }
801 58         192  
802             my($scales) = $in->copy->abs->maximum; # elementwise by rows
803 58 50       375  
804 58         514 if(($scales==0)->sum) {
805             return undef;
806             }
807              
808 0         0 # Some holding tanks
809 58 50       654 my($tmprow) = $out->slice('(0)')->zeroes;
810 58         417 $tmprow = $tmprow->double if $tmprow->type < double;
811             my($tmpval) = $tmprow->slice('(0)')->sever;
812 58         204  
813 58         147 my($col,$row);
814 58         142 for $col(0..$n1) {
815 125 100       286 for $row(1..$n1) {
816 160 100       450 my($klim) = $row<$col ? $row : $col;
817 160         501 if($klim > 0) {
818 92         137 $klim--;
819 92         1522 my($el) = $out->index2d($col,$row);
820             $el -= ( $out->slice("($col),0:$klim") *
821             $out->slice("0:$klim,($row)") )->sumover;
822             }
823              
824             }
825              
826             # Figure a_ij, with pivoting
827 92 100       423  
828             if($col < $n1) {
829 125         327 # Find the maximum value in the rest of the row
830 68         324 my $sl = $out->slice("($col),$col:$n1");
831 68         248 my $wh = $sl->abs->maximum_ind;
832             my $big = $sl->index($wh)->sever;
833              
834             # Permute if necessary to make the diagonal the maximum
835             # if($wh != 0)
836 68         1673 { # Permute rows to place maximum element on diagonal.
  68         130  
837             my $whc = $wh+$col;
838 68         280  
839 68         910 my $sl1 = $out->mv(1,0)->index($whc->slice("*$n"));
840 68         495 my $sl2 = $out->slice(":,($col)");
  68         246  
  68         201  
841             $tmprow .= $sl1; $sl1 .= $sl2; $sl2 .= $tmprow;
842 68         152  
843 68         918 $sl1 = $permute->index($whc);
844 68         816 $sl2 = $permute->index($col);
  68         173  
  68         153  
845             $tmpval .= $sl1; $sl1 .= $sl2; $sl2 .= $tmpval;
846 68         161  
847             $parity->where($wh>0) *= -1.0;
848             }
849              
850 68         307 # LAPACK cgetrf does not try fix singularity so nor do we, even though NR does
851 68 100       546 my $notbig = $big->where(abs($big) < $TINY);
852             return if !$notbig->isempty;
853              
854 68         372 # Divide by the diagonal element (which is now the largest element)
855 67         127 my $tout;
856             ($tout = $out->slice("($col),".($col+1).":$n1")) /= $big->slice('*1');
857             } # end of pivoting part
858             } # end of column loop
859 67 50       427  
860             wantarray ? ($out,$permute,$parity) : $out;
861             }
862              
863             EOD
864              
865             ######################################################################
866             pp_add_exported('','lu_decomp2');
867             pp_addpm(<<'EOD');
868              
869             =head2 lu_decomp2
870              
871             =for sig
872              
873             Signature: (a(m,m); [o]lu(m,m))
874              
875             =for ref
876              
877             LU decompose a matrix, with no row permutation
878              
879             =for usage
880              
881             ($lu, $perm, $parity) = lu_decomp2($x);
882              
883             $lu = lu_decomp2($x,$perm,$parity); # or
884             $lu = lu_decomp2($x); # $perm and $parity are optional
885              
886             lu_decomp($x->inplace,$perm,$parity); # or
887             lu_decomp($x->inplace); # $perm and $parity are optional
888              
889             =for description
890              
891             C works just like L, but it does B
892             pivoting at all. For compatibility with L, it
893             will give you a permutation list and a parity scalar if you ask
894             for them -- but they are always trivial.
895              
896             Because C does not pivot, it is numerically B --
897             that means it is less precise than L, particularly for
898             large or near-singular matrices. There are also specific types of
899             non-singular matrices that confuse it (e.g. ([0,-1,0],[1,0,0],[0,0,1]),
900             which is a 90 degree rotation matrix but which confuses C).
901              
902             On the other hand, if you want to invert rapidly a few hundred thousand
903             small matrices and don't mind missing one or two, it could be the ticket.
904             It can be up to 60% faster at the expense of possible failure of the
905             decomposition for some of the input matrices.
906              
907             The output is a single matrix that contains the LU decomposition of C<$a>;
908             you can even do it in-place, thereby destroying C<$a>, if you want. See
909             L for more information about LU decomposition.
910              
911             C is ported from I into PDL.
912              
913             =cut
914              
915             *PDL::lu_decomp2 = \&lu_decomp2;
916              
917 57     0 1 529 sub lu_decomp2 {
918 0         0 my($in) = shift;
919 0         0 my($perm) = shift;
920             my($par) = shift;
921 0         0  
922             my($sing_ok) = shift;
923 0         0  
924             my $TINY = 1e-30;
925 0 0 0     0  
      0        
926             barf("lu_decomp2 requires a square (2D) PDL\n")
927             if(!UNIVERSAL::isa($in,'PDL') ||
928             $in->ndims < 2 ||
929             $in->dim(0) != $in->dim(1));
930 0         0  
931 0         0 my($n) = $in->dim(0);
  0         0  
932             my($n1) = $n; $n1--;
933 0         0  
934 0 0       0 my($inplace) = $in->is_inplace;
935             my($out) = ($inplace) ? $in : $in->copy;
936 0 0       0  
937 0 0 0     0  
      0        
938             if(defined $perm) {
939             barf('lu_decomp2: permutation vector must match the matrix')
940             if(!UNIVERSAL::isa($perm,'PDL') ||
941 0         0 $perm->ndims != 1 ||
942             $perm->dim(0) != $out->dim(0));
943 0         0 $perm .= PDL->xvals($in->dim(0));
944             } else {
945             $perm = PDL->xvals($in->dim(0));
946 0 0       0 }
947 0 0 0     0  
948             if(defined $par) {
949             barf('lu_decomp: parity must be a scalar PDL')
950 0         0 if(!UNIVERSAL::isa($par,'PDL') ||
951             $par->nelem != 1);
952 0         0 $par .= 1.0;
953             } else {
954             $par = pdl(1.0);
955 0         0 }
956              
957 0         0 my $diagonal = $out->diagonal(0,1);
958 0         0  
959 0         0 my($col,$row);
960 0 0       0 for $col(0..$n1) {
961 0 0       0 for $row(1..$n1) {
962 0         0 my($klim) = $row<$col ? $row : $col;
963 0         0 if($klim > 0) {
964             $klim--;
965 0         0 my($el) = $out->index2d($col,$row);
966              
967             $el -= ( $out->slice("($col),0:$klim") *
968             $out->slice("0:$klim,($row)") )->sumover;
969             }
970              
971             }
972 0 0       0  
973             # Figure a_ij, with no pivoting
974 0         0 if($col < $n1) {
975             # Divide the rest of the column by the diagonal element
976             $out->slice("($col),".($col+1).":$n1") /= $diagonal->index($col)->dummy(0,$n1-$col);
977             }
978              
979 0 0       0 } # end of column loop
980              
981             wantarray ? ($out,$perm,$par) : $out;
982             }
983              
984             EOD
985              
986             ######################################################################
987             pp_add_exported('','lu_backsub');
988             pp_addpm(<<'EOD');
989              
990             =head2 lu_backsub
991              
992             =for sig
993              
994             Signature: (lu(m,m); perm(m); b(m))
995              
996             =for ref
997              
998             Solve A x = B for matrix A, by back substitution into A's LU decomposition.
999              
1000             =for usage
1001              
1002             ($lu,$perm,$par) = lu_decomp($A);
1003              
1004             $x = lu_backsub($lu,$perm,$par,$A); # or
1005             $x = lu_backsub($lu,$perm,$B); # $par is not required for lu_backsub
1006              
1007             lu_backsub($lu,$perm,$B->inplace); # modify $B in-place
1008              
1009             $x = lu_backsub(lu_decomp($A),$B); # (ignores parity value from lu_decomp)
1010              
1011             # starting from square matrix A and columns matrix B, with mathematically
1012             # correct dimensions
1013             $A = identity(4) + ones(4, 4);
1014             $A->slice('2,0') .= 0; # break symmetry to see if need transpose
1015             $B = sequence(2, 4);
1016             # all these functions take B as rows, interpret as though notional columns
1017             # mathematically confusing but can't change as back-compat and also
1018             # familiar to Fortran users, so just transpose inputs and outputs
1019              
1020             # using lu_backsub
1021             ($lu,$perm,$par) = lu_decomp($A);
1022             $x = lu_backsub($lu,$perm,$par, $B->transpose)->transpose;
1023              
1024             # using simq
1025             # remove all active dims
1026             @A_dims = $A->dims; @B_dims = $B->transpose->dims;
1027             splice @A_dims, 0, 2; splice @B_dims, 0, 1;
1028             @broadcast = PDL::Core::dims_filled(\@A_dims, \@B_dims);
1029             # simq modifies A, so need 1 copy per broadcast else non-first run has wrong A
1030             ($x) = simq($A->dupN(1,1,map +($A_dims[$_]//1)==1?$broadcast[$_]:1, 0..$#broadcast)->copy, $B->transpose, 0);
1031             $x = $x->inplace->transpose;
1032              
1033             # or with PDL::LinearAlgebra wrappers of LAPACK
1034             $x = msolve($A, $B);
1035              
1036             # or with LAPACK
1037             use PDL::LinearAlgebra::Real;
1038             getrf($lu=$A->copy, $ipiv=null, $info=null);
1039             getrs($lu, 1, $x=$B->transpose->copy, $ipiv, $info=null); # again, need transpose
1040             $x=$x->inplace->transpose;
1041              
1042             # or with GSL
1043             use PDL::GSL::LINALG;
1044             LU_decomp(my $lu=$A->copy, my $p=null, my $signum=null);
1045             # $B and $x, first dim is because GSL treats as vector, higher dims broadcast
1046             # so we transpose in and back out
1047             LU_solve($lu, $p, $B->transpose, my $x=null);
1048             $x=$x->inplace->transpose;
1049              
1050             # proof of the pudding is in the eating:
1051             print $A x $x;
1052              
1053             =for description
1054              
1055             Given the LU decomposition of a square matrix (from L),
1056             C does back substitution into the matrix to solve
1057             C for given vector C. It is separated from the
1058             C method so that you can call the cheap C
1059             multiple times and not have to do the expensive LU decomposition
1060             more than once.
1061              
1062             C acts on single vectors and broadcasts in the usual
1063             way, which means that it treats C<$y> as the I
1064             of the input. If you want to process a matrix, you must
1065             hand in the I of the matrix, and then transpose
1066             the output when you get it back. that is because pdls are
1067             indexed by (col,row), and matrices are (row,column) by
1068             convention, so a 1-D pdl corresponds to a row vector, not a
1069             column vector.
1070              
1071             If C<$lu> is dense and you have more than a few points to
1072             solve for, it is probably cheaper to find C with
1073             L, and just multiply C.) in fact,
1074             L works by calling C with the identity
1075             matrix.
1076              
1077             C is ported from section 2.3 of I.
1078             It is written in PDL but should probably be implemented in C.
1079              
1080             =cut
1081              
1082             *PDL::lu_backsub = \&lu_backsub;
1083              
1084 0     52 1 0 sub lu_backsub {
1085 52 50       269 my ($lu, $perm, $y, $par);
1086 52 100       128 print STDERR "lu_backsub: entering debug version...\n" if $PDL::debug;
    50          
1087 52         182 if(@_==3) {
1088             ($lu, $perm, $y) = @_;
1089 2         5 } elsif(@_==4) {
1090             ($lu, $perm, $par, $y) = @_;
1091             }
1092 50 50       117  
1093             barf("lu_backsub: LU decomposition is undef -- probably from a singular matrix.\n")
1094             unless defined($lu);
1095 52 50 33     137  
      33        
1096             barf("Usage: \$x = lu_backsub(\$lu,\$perm,\$y); all must be PDLs\n")
1097             unless(UNIVERSAL::isa($lu,'PDL') &&
1098             UNIVERSAL::isa($perm,'PDL') &&
1099             UNIVERSAL::isa($y,'PDL'));
1100 52         449  
1101 52         195 my $n = $y->dim(0);
  52         98  
1102             my $n1 = $n; $n1--;
1103              
1104             # Make sure broadcasting dimensions are compatible.
1105             # There are two possible sources of broadcast dims:
1106             #
1107             # (1) over multiple LU (i.e., $lu,$perm) instances
1108             # (2) over multiple B (i.e., $y) column instances
1109             #
1110             # The full dimensions of the function call looks like
1111             #
1112             # lu_backsub( lu(m,m,X), perm(m,X), b(m,Y) )
1113             #
1114             # where X is the list of extra LU dims and Y is
1115             # the list of extra B dims. We have several possible
1116             # cases:
1117             #
1118 52         114 # (1) Check that m dims are compatible
1119 52         227 my $ludims = pdl($lu->dims);
1120 52         288 my $permdims = pdl($perm->dims);
1121             my $bdims = pdl($y->dims);
1122 52 50       228  
1123             print STDERR "lu_backsub: called with args: \$lu$ludims, \$perm$permdims, \$y$bdims\n" if $PDL::debug;
1124 52         139  
1125 52 50 33     187 my $m = $ludims->slice("(0)"); # this is the sig dimension
      33        
      33        
1126             unless ( ($ludims->slice(0) == $m) and ($ludims->slice(1) == $m) and
1127 52         134 ($permdims->slice(0) == $m) and ($bdims->slice(0) == $m)) {
1128             barf "lu_backsub: mismatched sig dimensions";
1129             }
1130 0         0  
1131 52         954 my $lunumthr = $ludims->dim(0)-2;
1132 52         161 my $permnumthr = $permdims->dim(0)-1;
1133 52 50 33     151 my $bnumthr = $bdims->dim(0)-1;
1134 52         306 unless ( ($lunumthr == $permnumthr) and ($ludims->slice("1:-1") == $permdims)->all ) {
1135             barf "lu_backsub: \$lu and \$perm broadcast dims not equal! \n";
1136             }
1137              
1138 0 100 66     0 # (2) If X == Y then default broadcasting is ok
1139 52 50       717 if ( ($bnumthr==$permnumthr) and ($bdims==$permdims)->all) {
1140 4         15 print STDERR "lu_backsub: have explicit broadcast dims, goto BROADCAST_OK\n" if $PDL::debug;
1141             goto BROADCAST_OK;
1142             }
1143              
1144             # (3) If X == (x,Y) then add x dummy to lu,perm
1145              
1146             # (4) If ndims(X) > ndims(Y) then must have #3
1147              
1148             # (5) If ndims(X) < ndims(Y) then foreach
1149             # non-trivial leading dim in X (x0,x1,..)
1150             # insert dummy (x0,x1) into lu and perm
1151              
1152             # This means that broadcasting occurs over all
1153             # leading non-trivial (not length 1) dims of
1154             # B unless all the broadcast dims are explicitly
1155             # matched to the LU dims.
1156              
1157             BROADCAST_OK:
1158              
1159 4         26 # Permute the vector and make a copy if necessary.
1160 52 100       311 my $out = $y->dummy(1,$y->dim(0))->index($perm->dummy(1));
1161 52 50       1992 $out = $out->sever if !$y->is_inplace;
1162             print STDERR "lu_backsub: starting with \$out" . pdl($out->dims) . "\n" if $PDL::debug;
1163              
1164             # Make sure broadcasting over lu happens OK...
1165 52 50       161  
1166 52 0       297 if($out->ndims < $lu->ndims-1) {
1167 0         0 print STDERR "lu_backsub: adjusting dims for \$out" . pdl($out->dims) . "\n" if $PDL::debug;
1168 0         0 do {
1169             $out = $out->dummy(-1,$lu->dim($out->ndims+1));
1170             } while($out->ndims < $lu->ndims-1);
1171             }
1172              
1173 0         0 ## Do forward substitution into L
1174             my $row; my $r1;
1175 52         110  
1176 52         150 for $row(1..$n1) {
1177 57         109 $r1 = $row-1;
1178             $out->index($row) -= ($lu->slice("0:$r1,$row") *
1179             $out->slice("0:$r1")
1180             )->sumover;
1181             }
1182              
1183 57         877 ## Do backward substitution into U, and normalize by the diagonal
1184 52         178 my $ludiag = $lu->diagonal(0,1);
1185             $out->index($n1) /= $ludiag->index($n1)->dummy(0); # TODO: check broadcasting
1186 52         857  
1187 52         681 for ($row=$n1; $row>0; $row--) {
1188 57         117 $r1 = $row-1;
1189             $out->index($r1) -= ($lu->slice("$row:$n1,$r1") * # TODO: check broadcast dims
1190             $out->slice("$row:$n1")
1191 57         810 )->sumover;
1192             $out->index($r1) /= $ludiag->index($r1)->dummy(0); # TODO: check broadcast dims
1193             }
1194 57 100       1570  
1195 52 100       195 if ($y->is_inplace) {
1196 3         11 $y->setdims([$out->dims]) if !PDL::all($y->shape == $out->shape); # assgn needs same shape
1197             $y .= $out;
1198 3         38 }
1199             $out;
1200             }
1201             EOD
1202              
1203             # XXX Destroys a!!!
1204             # To use the new a again, must store both a and ips.
1205             pp_def("simq",
1206             HandleBad => 0,
1207             Pars => '[io,phys]a(n,n); [phys]b(n); [o,phys]x(n); int [o,phys]ips(n)',
1208             OtherPars => 'int flag;',
1209             GenericTypes => ['D'],
1210             Code => '
1211             extern int simq( double *A, double *B, double *X,
1212             int n, int flag, int *IPS );
1213             simq($P(a), $P(b), $P(x), $SIZE(n), $COMP(flag), $P(ips));
1214             ',
1215             Doc => '
1216             =for ref
1217              
1218             Solution of simultaneous linear equations, C.
1219             B.
1220              
1221             =for usage
1222              
1223             # remove all active dims
1224             @A_dims = $A->dims; @B_dims = $B->transpose->dims;
1225             splice @A_dims, 0, 2; splice @B_dims, 0, 1;
1226             @broadcast = PDL::Core::dims_filled(\@A_dims, \@B_dims);
1227             # simq modifies A, so need 1 copy per broadcast else non-first run has wrong A
1228             ($x) = simq($A->dupN(1,1,map +($A_dims[$_]//1)==1?$broadcast[$_]:1, 0..$#broadcast)->copy, $B->transpose, 0);
1229             $x = $x->inplace->transpose;
1230              
1231             C<$a> is an C matrix (i.e., a vector of length C), stored row-wise:
1232             that is, C, where C.
1233              
1234             While this is the transpose of the normal column-wise storage, this
1235             corresponds to normal PDL usage. The contents of matrix a may be
1236             altered (but may be required for subsequent calls with flag = -1).
1237              
1238             C<$y>, C<$x>, C<$ips> are vectors of length C.
1239              
1240             Set C to solve.
1241             Set C to do a new back substitution for
1242             different C<$y> vector using the same a matrix previously reduced when
1243             C (the C<$ips> vector generated in the previous solution is also
1244             required).
1245              
1246             For this function to work well with broadcasting, it will need the LU
1247             decomposition part split out, so that for solving C only C
1248             would be written to.
1249              
1250             See also L, which does the same thing with a slightly
1251             less opaque interface.
1252              
1253             =cut
1254              
1255              
1256             ');
1257              
1258             ######################################################################
1259             ### squaretotri
1260             ###
1261             # this doesn't need to be changed to support bad values
1262             # I could put 'HandleBad => 1', but it would just cause an
1263             # unnecessary increase (admittedly small) in the amount of
1264             # code
1265             #
1266             pp_def("squaretotri",
1267             Pars => 'a(n,n); [o]b(m=CALC(($SIZE(n) * ($SIZE(n)+1))/2))',
1268             GenericTypes => [ppdefs_all],
1269             Code => '
1270             register PDL_Indx mna=0, nb=0;
1271             loop(m) %{
1272             $b() = $a(n0 => mna, n1 => nb);
1273             mna++; if(mna > nb) {mna = 0; nb ++;}
1274             %}
1275             ',
1276             Doc => '=for ref
1277              
1278             Convert a lower-triangular square matrix to triangular vector storage.
1279             Ignores upper half of input.
1280             ',
1281             );
1282              
1283             pp_def("tritosquare",
1284             Pars => 'a(m); [o]b(n,n)',
1285             GenericTypes => [ppdefs_all],
1286             RedoDimsCode => $TRIANGULAR_REDODIMS,
1287             Code => '
1288             register PDL_Indx mna=0, nb=0;
1289             loop(m) %{
1290             $b(n0 => mna, n1 => nb) = $a();
1291             mna++; if(mna > nb) {mna = 0; nb ++;}
1292             %}
1293             ',
1294             Doc => '=for ref
1295              
1296             Convert a triangular vector to lower-triangular square matrix storage.
1297             Does not touch upper half of output.
1298             ',
1299             );
1300              
1301             pp_def('tricpy',
1302             Pars => 'A(m,n);[o] C(m,n)',
1303             OtherPars => 'int uplo',
1304             OtherParsDefaults => {uplo => 0},
1305             ArgOrder => [qw(A uplo C)],
1306             GenericTypes => [ppdefs_all()],
1307             Code => '
1308             if ($COMP(uplo)) { /* lower */
1309             broadcastloop %{ loop(n,m=:n+1) %{ $C() = $A(); %} %}
1310             } else {
1311             broadcastloop %{ loop(n,m=n) %{ $C() = $A(); %} %}
1312             }
1313             ',
1314             Doc => <<'EOT',
1315             =for ref
1316              
1317             Copy triangular part to another matrix. If uplo == 0 copy upper triangular
1318             part.
1319              
1320             Originally by Grégory Vanuxem.
1321             EOT
1322             );
1323              
1324             pp_def('mstack',
1325             DefaultFlow => 1,
1326             TwoWay => 1,
1327             Pars => 'x(n,m);y(n,p);[o]out(n,q=CALC($SIZE(m)+$SIZE(p)));',
1328             GenericTypes => [ppdefs_all()],
1329             Code => '
1330             loop (m,n) %{ $out(q=>m) = $x(); %}
1331             loop (q=$SIZE(m),n) %{ $out() = $y(p=>q-$SIZE(m)); %}
1332             ',
1333             BackCode => '
1334             loop (m,n) %{ $x() = $out(q=>m); %}
1335             loop (q=$SIZE(m),n) %{ $y(p=>q-$SIZE(m)) = $out(); %}
1336             ',
1337             Doc => <
1338             =for ref
1339              
1340             Combine two 2D ndarrays into a single ndarray, along the second
1341             ("vertical") dim.
1342             This routine does backward and forward dataflow automatically.
1343              
1344             Originally by Grégory Vanuxem.
1345             EOT
1346             );
1347              
1348             pp_def('augment',
1349             DefaultFlow => 1,
1350             TwoWay => 1,
1351             Pars => 'x(n); y(p);[o]out(q=CALC($SIZE(n)+$SIZE(p)))',
1352             GenericTypes => [ppdefs_all()],
1353             Code => '
1354             loop (q=:$SIZE(n)) %{ $out() = $x(n=>q); %}
1355             loop (q=$SIZE(n)) %{ $out() = $y(p=>q-$SIZE(n)); %}
1356             ',
1357             BackCode => '
1358             loop (q=:$SIZE(n)) %{ $x(n=>q) = $out(); %}
1359             loop (q=$SIZE(n)) %{ $y(p=>q-$SIZE(n)) = $out(); %}
1360             ',
1361             Doc => <
1362             =for ref
1363              
1364             Combine two ndarrays into a single ndarray along the 0-th ("horizontal") dim.
1365             This routine does backward and forward dataflow automatically.
1366              
1367             Originally by Grégory Vanuxem.
1368             EOT
1369             );
1370              
1371             pp_addpm({At=>'Bot'},<<'EOD');
1372              
1373             =head1 AUTHOR
1374              
1375             Copyright (C) 2002 Craig DeForest (deforest@boulder.swri.edu),
1376             R.J.R. Williams (rjrw@ast.leeds.ac.uk), Karl Glazebrook
1377             (kgb@aaoepp.aao.gov.au). There is no warranty. You are allowed to
1378             redistribute and/or modify this work under the same conditions as PDL
1379             itself. If this file is separated from the PDL distribution, then the
1380             PDL copyright notice should be included in this file.
1381              
1382             =cut
1383              
1384             EOD
1385              
1386             pp_done();