File Coverage

blib/lib/PDLA/MatrixOps.pm
Criterion Covered Total %
statement 202 269 75.0
branch 62 134 46.2
condition 27 90 30.0
subroutine 15 18 83.3
pod 8 11 72.7
total 314 522 60.1


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