2 * Introduction to lapack::
3 * Functions and Variables for lapack::
6 @node Introduction to lapack, Functions and Variables for lapack, Package lapack, Package lapack
7 @section Introduction to lapack
9 @code{lapack} is a Common Lisp translation (via the program @code{f2cl}) of the Fortran library LAPACK,
10 as obtained from the SLATEC project.
12 @opencatbox{Categories:}
13 @category{Numerical methods}
14 @category{Share packages}
15 @category{Package lapack}
19 @node Functions and Variables for lapack, , Introduction to lapack, Package lapack
20 @section Functions and Variables for lapack
23 @deffn {Function} dgeev (@var{A})
24 @deffnx {Function} dgeev (@var{A}, @var{right_p}, @var{left_p})
26 Computes the eigenvalues and, optionally, the eigenvectors of a matrix @var{A}.
27 All elements of @var{A} must be integer or floating point numbers.
28 @var{A} must be square (same number of rows and columns).
29 @var{A} might or might not be symmetric.
31 To make use of this function, you must load the LaPack package via
32 @code{load("lapack")}.
34 @code{dgeev(@var{A})} computes only the eigenvalues of @var{A}.
35 @code{dgeev(@var{A}, @var{right_p}, @var{left_p})} computes the eigenvalues of @var{A}
36 and the right eigenvectors when @math{@var{right_p} = @code{true}}
37 and the left eigenvectors when @math{@var{left_p} = @code{true}}.
39 A list of three items is returned.
40 The first item is a list of the eigenvalues.
41 The second item is @code{false} or the matrix of right eigenvectors.
42 The third item is @code{false} or the matrix of left eigenvectors.
45 m4_math(<<<v_j>>>,<<<v_j>>>)
46 (the @math{j}-th column of the right eigenvector matrix) satisfies
49 <<<\mathbf{A} v_j = \lambda_j v_j>>>,
50 <<<@math{A . v_j = lambda_j . v_j}>>>,
51 <<<{\bf A} v_j = \lambda_j v_j>>>
55 m4_math(<<<\lambda_j>>>,<<<@math{lambda_j}>>>)
56 is the corresponding eigenvalue.
58 m4_math(<<<u_j>>>,<<<u_j>>>)
59 (the @math{j}-th column of the left eigenvector matrix) satisfies
62 <<<u_j^\mathbf{H} \mathbf{A} = \lambda_j u_j^\mathbf{H}>>>,
63 <<<@math{u_j^H . A = lambda_j . u_j^H}>>>,
64 <<<u_j^{\bf H} {\bf A} = \lambda_j u_j^{\bf H}>>>
68 m4_math(<<<u_j^\mathbf{H}>>>,<<<u_j^H>>>, <<<u_j^{\bf H}>>>)
69 denotes the conjugate transpose of
70 m4_mathdot(<<<u_j>>>,<<<u_j>>>)
71 For a Maxima function to compute the conjugate transpose, @pxref{ctranspose}.
74 The computed eigenvectors are normalized to have Euclidean norm
75 equal to 1, and largest component has imaginary part equal to zero.
82 @c M : matrix ([9.5, 1.75], [3.25, 10.45]);
84 @c [L, v, u] : dgeev (M, true, true);
85 @c D : apply (diag_matrix, L);
87 @c transpose (u) . M - D . transpose (u);
90 (%i1) load ("lapack")$
92 (%i2) fpprintprec : 6;
96 (%i3) M : matrix ([9.5, 1.75], [3.25, 10.45]);
103 (%o4) [[7.54331, 12.4067], false, false]
106 (%i5) [L, v, u] : dgeev (M, true, true);
107 [ - 0.666642 - 0.515792 ]
108 (%o5) [[7.54331, 12.4067], [ ],
109 [ 0.745378 - 0.856714 ]
110 [ - 0.856714 - 0.745378 ]
112 [ 0.515792 - 0.666642 ]
115 (%i6) D : apply (diag_matrix, L);
122 [ 0.0 - 8.88178e-16 ]
124 [ - 8.88178e-16 0.0 ]
127 (%i8) transpose (u) . M - D . transpose (u);
128 [ 0.0 - 4.44089e-16 ]
134 @opencatbox{Categories:}
135 @category{Package lapack}
141 @deffn {Function} dgeqrf (@var{A})
143 Computes the QR decomposition of the matrix @var{A}.
144 All elements of @var{A} must be integer or floating point numbers.
145 @var{A} may or may not have the same number of rows and columns.
147 To make use of this function, you must load the LaPack package via
148 @code{load("lapack")}.
150 The real square matrix
151 m4_math(<<<\mathbf{A}>>>, <<<A>>>, <<<{\bf A}>>>)
155 <<<\mathbf{A} = \mathbf{Q}\mathbf{R}>>>,
157 <<<{\bf A} = {\bf Q} {\bf R}>>>)
160 m4_math(<<<{\bf Q}>>>, <<<Q>>>)
161 is a square orthogonal matrix with the same number of rows as
162 m4_math(<<<\mathbf{A}>>>, <<<A>>>, <<<{\bf A}>>>)
164 m4_math(<<<{\bf R}>>>, <<<R>>>)
165 is an upper triangular matrix and is the same size as
166 m4_mathdot(<<<{\bf A}>>>, A)
168 A list of two items is returned.
169 The first item is the matrix
170 m4_mathdot(<<<{\bf Q}>>>, Q)
171 The second item is the matrix
172 m4_mathcomma(<<<{\bf R}>>>, R)
173 The product @math{Q . R}, where "." is the noncommutative multiplication operator,
174 is equal to @var{A} (ignoring floating point round-off errors).
179 @c M : matrix ([1, -3.2, 8], [-11, 2.7, 5.9]);
180 @c [q, r] : dgeqrf (M);
185 (%i1) load ("lapack")$
187 (%i2) fpprintprec : 6;
191 (%i3) M : matrix ([1, -3.2, 8], [-11, 2.7, 5.9]);
197 (%i4) [q, r] : dgeqrf (M);
198 [ - 0.0905357 0.995893 ]
200 [ 0.995893 0.0905357 ]
201 [ - 11.0454 2.97863 5.15148 ]
203 [ 0.0 - 2.94241 8.50131 ]
207 [ - 7.77156e-16 1.77636e-15 - 8.88178e-16 ]
209 [ 0.0 - 1.33227e-15 8.88178e-16 ]
212 (%i6) mat_norm (%, 1);
217 @opencatbox{Categories:}
218 @category{Package lapack}
224 @deffn {Function} dgesv (@var{A}, @var{b})
226 Computes the solution @var{x} of the linear equation
227 m4_mathcomma(<<<{\bf A} x = b>>>, <<<@var{A} @var{x} = @var{b}>>>)
229 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
230 is a square matrix, and @math{b} is a matrix of the same number of rows
232 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
233 and any number of columns.
234 The return value @math{x} is the same size as @math{b}.
236 To make use of this function, you must load the LaPack package via
237 @code{load("lapack")}.
239 The elements of @var{A} and @var{b} must evaluate to real floating point numbers via @code{float};
240 thus elements may be any numeric type, symbolic numerical constants, or expressions which evaluate to floats.
241 The elements of @var{x} are always floating point numbers.
242 All arithmetic is carried out as floating point operations.
244 @code{dgesv} computes the solution via the LU decomposition of @var{A}.
248 @code{dgesv} computes the solution of the linear equation @math{@var{A} @var{x} = @var{b}}.
251 @c A : matrix ([1, -2.5], [0.375, 5]);
252 @c b : matrix ([1.75], [-0.625]);
254 @c dlange (inf_norm, b - A . x);
258 (%i1) A : matrix ([1, -2.5], [0.375, 5]);
264 (%i2) b : matrix ([1.75], [-0.625]);
270 (%i3) x : dgesv (A, b);
272 (%o3) dgesv([ ], [ ])
273 [ 0.375 5 ] [ - 0.625 ]
276 (%i4) dlange (inf_norm, b - A . x);
278 (%o4) dlange(inf_norm, [ ]
280 [ 1 - 2.5 ] [ 1 - 2.5 ] [ 1.75 ]
281 - [ ] . dgesv([ ], [ ]))
282 [ 0.375 5 ] [ 0.375 5 ] [ - 0.625 ]
286 @var{b} is a matrix with the same number of rows as @var{A} and any number of columns.
287 @var{x} is the same size as @var{b}.
290 @c A : matrix ([1, -0.15], [1.82, 2]);
291 @c b : matrix ([3.7, 1, 8], [-2.3, 5, -3.9]);
293 @c dlange (inf_norm, b - A . x);
297 (%i1) A : matrix ([1, -0.15], [1.82, 2]);
303 (%i2) b : matrix ([3.7, 1, 8], [-2.3, 5, -3.9]);
309 (%i3) x : dgesv (A, b);
310 [ 1 - 0.15 ] [ 3.7 1 8 ]
311 (%o3) dgesv([ ], [ ])
312 [ 1.82 2 ] [ - 2.3 5 - 3.9 ]
315 (%i4) dlange (inf_norm, b - A . x);
317 (%o4) dlange(inf_norm, [ ]
319 [ 1 - 0.15 ] [ 1 - 0.15 ]
321 [ 1.82 2 ] [ 1.82 2 ]
328 The elements of @var{A} and @var{b} must evaluate to real floating point numbers.
331 @c A : matrix ([5, -%pi], [1b0, 11/17]);
332 @c b : matrix ([%e], [sin(1)]);
334 @c dlange (inf_norm, b - A . x);
338 (%i1) A : matrix ([5, -%pi], [1b0, 11/17]);
346 (%i2) b : matrix ([%e], [sin(1)]);
352 (%i3) x : dgesv (A, b);
355 (%o3) dgesv([ 11 ], [ ])
356 [ 1.0b0 -- ] [ sin(1) ]
360 (%i4) dlange (inf_norm, b - A . x);
362 (%o4) dlange(inf_norm, [ ]
364 [ 5 - %pi ] [ 5 - %pi ]
366 - [ 11 ] . dgesv([ 11 ], [ ]))
367 [ 1.0b0 -- ] [ 1.0b0 -- ] [ sin(1) ]
372 @opencatbox{Categories:}
373 @category{Package lapack}
374 @category{Linear equations}
380 @deffn {Function} dgesvd (@var{A})
381 @deffnx {Function} dgesvd (@var{A}, @var{left_p}, @var{right_p})
383 Computes the singular value decomposition (SVD) of a matrix @var{A},
384 comprising the singular values and, optionally, the left and right singular vectors.
385 All elements of @var{A} must be integer or floating point numbers.
386 @var{A} might or might not be square (same number of rows and columns).
388 To make use of this function, you must load the LaPack package via
389 @code{load("lapack")}.
391 Let @math{m} be the number of rows, and @math{n} the number of columns of @var{A}.
392 The singular value decomposition of
393 m4_math(<<<\mathbf{A}>>>,<<<@math{A}>>>, <<<{\bf A}>>>)
394 comprises three matrices,
395 m4_mathcomma(<<<\mathbf{U}>>>,<<<U>>>,<<<{\bf U}>>>)
396 m4_mathcomma(<<<\mathbf{\Sigma}>>>,<<<@math{Sigma}>>>, <<<{\bf \Sigma}>>>)
398 m4_mathcomma(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
401 @c this code breaks texi2pdf: @math{@var{A} = @var{U} . @var{Sigma} . @var{V^T}}
402 @c @math{@var{A} = @var{U} . @var{Sigma} . @var{V}^T}
405 <<<\mathbf{A} = \mathbf{U} \mathbf{\Sigma} \mathbf{V}^T>>>,
406 <<<@math{{A} = {U} . {Sigma} . {V}^T}>>>,
407 <<<{\bf A} = {\bf U} {\bf \Sigma} {\bf V}^T>>>)
411 m4_math(<<<\mathbf{U}>>>, <<<@math{U}>>>, <<<{\bf U}>>>)
412 is an @math{m}-by-@math{m} unitary matrix,
413 m4_math(<<<\mathbf{\Sigma}>>>, <<<@math{Sigma}>>>, <<<{\bf\Sigma}>>>)
414 is an @math{m}-by-@math{n} diagonal matrix,
416 m4_math(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
417 is an @math{n}-by-@math{n} unitary matrix.
420 m4_math(<<<\mathbf{\sigma}_i>>>, <<<@math{sigma[i]}>>>, <<<{\bf \sigma}_i>>>)
421 be a diagonal element of
422 m4_mathcomma(<<<\mathbf{\Sigma}>>>, <<<@math{Sigma}>>>, <<<{\bf \Sigma}>>>)
424 m4_mathdot(<<<\mathbf{\Sigma}_{ii} = \sigma_i>>>,
425 <<<@math{@var{Sigma}[i, i] = @var{sigma}[i]}>>>,
426 <<<{\bf \Sigma}_{ii} = \sigma_i>>>)
428 m4_math(<<<\sigma_i>>>, <<<@math{sigma[i]}>>>)
429 are the so-called singular values of
430 m4_mathpunc(<<<;>>>, <<<\mathbf{A}>>>, <<<A>>>, <<<{\bf A}>>>)
431 these are real and nonnegative, and returned in descending order.
433 m4_math(<<<\min(m, n)>>>, <<<min(m, n)>>>)
435 m4_math(<<<\mathbf{U}>>>, <<<@math{U}>>>, <<<{\bf U}>>>)
437 m4_math(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
438 are the left and right singular vectors of
439 m4_mathdot(<<<\mathbf{A}>>>, <<<@math{A}>>>, <<<{\bf A}>>>)
440 Note that @code{dgesvd} returns the transpose of
441 m4_mathcomma(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
443 m4_math(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
446 @code{dgesvd(@var{A})} computes only the singular values of @var{A}.
447 @code{dgesvd(@var{A}, @var{left_p}, @var{right_p})} computes the singular values of @var{A}
448 and the left singular vectors when @math{@var{left_p} = @code{true}}
449 and the right singular vectors when @math{@var{right_p} = @code{true}}.
451 A list of three items is returned.
452 The first item is a list of the singular values.
453 The second item is @code{false} or the matrix of left singular vectors.
454 The third item is @code{false} or the matrix of right singular vectors.
461 @c M: matrix([1, 2, 3], [3.5, 0.5, 8], [-1, 2, -3], [4, 9, 7]);
463 @c [sigma, U, VT] : dgesvd (M, true, true);
467 @c genmatrix(lambda ([i, j], if i=j then sigma[i] else 0),
469 @c U . Sigma . VT - M;
470 @c transpose (U) . U;
471 @c VT . transpose (VT);
474 (%i1) load ("lapack")$
476 (%i2) fpprintprec : 6;
480 (%i3) M: matrix([1, 2, 3], [3.5, 0.5, 8], [-1, 2, -3], [4, 9, 7]);
491 (%o4) [[14.4744, 6.38637, 0.452547], false, false]
494 (%i5) [sigma, U, VT] : dgesvd (M, true, true);
495 (%o5) [[14.4744, 6.38637, 0.452547],
496 [ - 0.256731 0.00816168 0.959029 - 0.119523 ]
498 [ - 0.526456 0.672116 - 0.206236 - 0.478091 ]
500 [ 0.107997 - 0.532278 - 0.0708315 - 0.83666 ]
502 [ - 0.803287 - 0.514659 - 0.180867 0.239046 ]
503 [ - 0.374486 - 0.538209 - 0.755044 ]
505 [ 0.130623 - 0.836799 0.5317 ]]
507 [ - 0.917986 0.100488 0.383672 ]
510 (%i6) m : length (U);
514 (%i7) n : length (VT);
519 genmatrix(lambda ([i, j], if i=j then sigma[i] else 0),
530 (%i9) U . Sigma . VT - M;
531 [ 1.11022e-15 0.0 1.77636e-15 ]
533 [ 1.33227e-15 1.66533e-15 0.0 ]
535 [ - 4.44089e-16 - 8.88178e-16 4.44089e-16 ]
537 [ 8.88178e-16 1.77636e-15 8.88178e-16 ]
540 (%i10) transpose (U) . U;
541 [ 1.0 5.55112e-17 2.498e-16 2.77556e-17 ]
543 [ 5.55112e-17 1.0 5.55112e-17 4.16334e-17 ]
545 [ 2.498e-16 5.55112e-17 1.0 - 2.08167e-16 ]
547 [ 2.77556e-17 4.16334e-17 - 2.08167e-16 1.0 ]
550 (%i11) VT . transpose (VT);
551 [ 1.0 0.0 - 5.55112e-17 ]
553 (%o11) [ 0.0 1.0 5.55112e-17 ]
555 [ - 5.55112e-17 5.55112e-17 1.0 ]
559 @opencatbox{Categories:}
560 @category{Package lapack}
567 @deffn {Function} dlange (@var{norm}, @var{A})
568 @deffnx {Function} zlange (@var{norm}, @var{A})
570 Computes a norm or norm-like function of the matrix @var{A}. If
571 @var{A} is a real matrix, use @code{dlange}. For a matrix with
572 complex elements, use @code{zlange}.
574 To make use of this function, you must load the LaPack package via
575 @code{load("lapack")}.
577 @code{norm} specifies the kind of norm to be computed:
581 m4_math(<<<\max(|{\bf A}_{ij}|)>>>, <<<max(abs(A(i, j)))>>>)
582 where @math{i} and @math{j} range over
583 the rows and columns, respectively, of
584 m4_mathdot(<<<{\bf A}>>>, <<<A>>>)
585 Note that this function is not a proper matrix norm.
588 m4_math(<<<L_1>>>, <<<L[1]>>>)
590 m4_mathcomma(<<<{\bf A}>>>, <<<A>>>)
591 that is, the maximum of the sum of the absolute value of elements in each column.
594 m4_math(<<<L_\infty>>>, <<<L[inf]>>>)
596 m4_mathcomma(<<<{\bf A}>>>, <<<A>>>)
597 that is, the maximum of the sum of the absolute value of elements in each row.
599 Compute the Frobenius norm of
600 m4_mathcomma(<<<{\bf A}>>>, <<<A>>>)
601 that is, the square root of the sum of squares of the matrix elements.
604 @opencatbox{Categories:}
605 @category{Package lapack}
611 @deffn {Function} dgemm (@var{A}, @var{B})
612 @deffnx {Function} dgemm (@var{A}, @var{B}, @var{options})
613 Compute the product of two matrices and optionally add the product to
616 In the simplest form, @code{dgemm(@var{A}, @var{B})} computes the
617 product of the two real matrices, @var{A} and @var{B}.
619 To make use of this function, you must load the LaPack package via
620 @code{load("lapack")}.
622 In the second form, @code{dgemm} computes
623 m4_math(<<<\alpha {\bf A} {\bf B} + \beta {\bf C}>>>,
624 <<<@var{alpha} * @var{A} * @var{B} + @var{beta} * @var{C}>>>)
626 m4_mathcomma(<<<{\bf A}>>>, <<<@var{A}>>>)
627 m4_mathcomma(<<<{\bf B}>>>, <<<@var{B}>>>)
629 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
630 are real matrices of the appropriate sizes and
631 m4_math(<<<\alpha>>>, <<<alpha>>>)
633 m4_math(<<<\beta>>>, <<<beta>>>)
634 are real numbers. Optionally,
635 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
637 m4_math(<<<{\bf B}>>>, <<<@var{B}>>>)
639 be transposed before computing the product. The extra parameters are
640 specified by optional keyword arguments: The keyword arguments are
641 optional and may be specified in any order. They all take the form
642 @code{key=val}. The keyword arguments are:
647 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
648 that should be added. The default is @code{false},
649 which means no matrix is added.
652 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
654 m4_math(<<<{\bf B}>>>, <<<@var{B}>>>)
655 is multiplied by this value. The
659 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
660 is given, this value multiplies
661 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
663 is added. The default value is 0, which implies that
664 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
667 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
668 is given. Hence, be sure to specify a non-zero
670 m4_mathdot(<<<\beta>>>, <<<@math{beta}>>>)
672 If @code{true}, the transpose of
673 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
675 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
676 for the product. The default is @code{false}.
678 If @code{true}, the transpose of
679 m4_math(<<<{\bf B}>>>, <<<@var{B}>>>)
681 m4_math(<<<{\bf B}>>>, <<<@var{B}>>>)
682 for the product. The default is @code{false}.
687 @c A : matrix([1,2,3],[4,5,6],[7,8,9]);
688 @c B : matrix([-1,-2,-3],[-4,-5,-6],[-7,-8,-9]);
689 @c C : matrix([3,2,1],[6,5,4],[9,8,7]);
692 @c dgemm(A,B,transpose_a=true);
694 @c dgemm(A,B,c=C,beta=1);
696 @c dgemm(A,B,c=C,beta=1, alpha=-1);
700 (%i1) load ("lapack")$
702 (%i2) A : matrix([1,2,3],[4,5,6],[7,8,9]);
710 (%i3) B : matrix([-1,-2,-3],[-4,-5,-6],[-7,-8,-9]);
713 (%o3) [ - 4 - 5 - 6 ]
718 (%i4) C : matrix([3,2,1],[6,5,4],[9,8,7]);
727 [ - 30.0 - 36.0 - 42.0 ]
729 (%o5) [ - 66.0 - 81.0 - 96.0 ]
731 [ - 102.0 - 126.0 - 150.0 ]
737 (%o6) [ - 66 - 81 - 96 ]
739 [ - 102 - 126 - 150 ]
742 (%i7) dgemm(A,B,transpose_a=true);
743 [ - 66.0 - 78.0 - 90.0 ]
745 (%o7) [ - 78.0 - 93.0 - 108.0 ]
747 [ - 90.0 - 108.0 - 126.0 ]
750 (%i8) transpose(A) . B;
753 (%o8) [ - 78 - 93 - 108 ]
758 (%i9) dgemm(A,B,c=C,beta=1);
759 [ - 27.0 - 34.0 - 41.0 ]
761 (%o9) [ - 60.0 - 76.0 - 92.0 ]
763 [ - 93.0 - 118.0 - 143.0 ]
769 (%o10) [ - 60 - 76 - 92 ]
774 (%i11) dgemm(A,B,c=C,beta=1, alpha=-1);
777 (%o11) [ 72.0 86.0 100.0 ]
779 [ 111.0 134.0 157.0 ]
790 @opencatbox{Categories:}
791 @category{Package lapack}
797 @deffn {Function} zgeev (@var{A})
798 @deffnx {Function} zgeev (@var{A}, @var{right_p}, @var{left_p})
800 Like @mrefcomma{dgeev} but the matrix
801 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
804 To make use of this function, you must load the LaPack package via
805 @code{load("lapack")}.
807 @opencatbox{Categories:}
808 @category{Package lapack}
814 @deffn {Function} zheev (@var{A})
815 @deffnx {Function} zheev (@var{A}, @var{eigvec_p})
817 Like @mrefcomma{dgeev} but the matrix
818 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
819 is assumed to be a square
820 complex Hermitian matrix. If @var{eigvec_p} is @code{true}, then the
821 eigenvectors of the matrix are also computed.
823 To make use of this function, you must load the LaPack package via
824 @code{load("lapack")}.
826 No check is made that the matrix
827 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
828 is, in fact, Hermitian.
830 A list of two items is returned, as in @code{dgeev}: a list of
831 eigenvalues, and @code{false} or the matrix of the eigenvectors.
836 @c [9.14 +%i*0.00 , -4.37 -%i*9.22 , -1.98 -%i*1.72 , -8.96 -%i*9.50],
837 @c [-4.37 +%i*9.22 , -3.35 +%i*0.00 , 2.25 -%i*9.51 , 2.57 +%i*2.40],
838 @c [-1.98 +%i*1.72 , 2.25 +%i*9.51 , -4.82 +%i*0.00 , -3.24 +%i*2.04],
839 @c [-8.96 +%i*9.50 , 2.57 -%i*2.40 , -3.24 -%i*2.04 , 8.44 +%i*0.00]);
846 (%i1) load("lapack")$
849 [9.14 +%i*0.00 , -4.37 -%i*9.22 , -1.98 -%i*1.72 , -8.96 -%i*9.50],
850 [-4.37 +%i*9.22 , -3.35 +%i*0.00 , 2.25 -%i*9.51 , 2.57 +%i*2.40],
851 [-1.98 +%i*1.72 , 2.25 +%i*9.51 , -4.82 +%i*0.00 , -3.24 +%i*2.04],
852 [-8.96 +%i*9.50 , 2.57 -%i*2.40 , -3.24 -%i*2.04 , 8.44 +%i*0.00]);
853 [ 9.14 ] [ - 9.22 %i - 4.37 ]
855 [ 9.22 %i - 4.37 ] [ - 3.35 ]
856 (%o2) Col 1 = [ ] Col 2 = [ ]
857 [ 1.72 %i - 1.98 ] [ 9.51 %i + 2.25 ]
859 [ 9.5 %i - 8.96 ] [ 2.57 - 2.4 %i ]
860 [ - 1.72 %i - 1.98 ] [ - 9.5 %i - 8.96 ]
862 [ 2.25 - 9.51 %i ] [ 2.4 %i + 2.57 ]
863 Col 3 = [ ] Col 4 = [ ]
864 [ - 4.82 ] [ 2.04 %i - 3.24 ]
866 [ - 2.04 %i - 3.24 ] [ 8.44 ]
870 (%o3) [[- 16.004746472094734, - 6.764970154793324,
871 6.6657114535070985, 25.51400517338097], false]
873 (%i4) E: zheev(M,true)$
876 (%o5) [- 16.004746472094737, - 6.764970154793325,
877 6.665711453507101, 25.514005173380962]
881 [ 0.26746505331727455 %i + 0.21754535866650165 ]
883 [ 0.002696730886619885 %i + 0.6968836773391712 ]
885 [ - 0.6082406376714117 %i - 0.012106142926979313 ]
887 [ 0.15930818580950368 ]
888 [ 0.26449374706674444 %i + 0.4773693349937472 ]
890 [ - 0.28523890360316206 %i - 0.14143627420116733 ]
892 [ 0.2654607680986639 %i + 0.44678181171841735 ]
894 [ 0.5750762708542709 ]
895 [ 0.28106497673059216 %i - 0.13352639282451817 ]
897 [ 0.28663101328695556 %i - 0.4536971347853274 ]
899 [ - 0.29336843237542953 %i - 0.49549724255410565 ]
901 [ 0.5325337537576771 ]
902 [ - 0.5737316575503476 %i - 0.39661467994277055 ]
904 [ 0.018265026190214573 %i + 0.35305577043870173 ]
906 [ 0.16737009000854253 %i + 0.01476684746229564 ]
908 [ 0.6002632636961784 ]
912 @opencatbox{Categories:}
913 @category{Package lapack}