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")$
91 (%i2) fpprintprec : 6;
93 (%i3) M : matrix ([9.5, 1.75], [3.25, 10.45]);
98 (%o4) [[7.54331, 12.4067], false, false]
99 (%i5) [L, v, u] : dgeev (M, true, true);
100 [ - .666642 - .515792 ]
101 (%o5) [[7.54331, 12.4067], [ ],
102 [ .745378 - .856714 ]
103 [ - .856714 - .745378 ]
105 [ .515792 - .666642 ]
106 (%i6) D : apply (diag_matrix, L);
111 [ 0.0 - 8.88178E-16 ]
113 [ - 8.88178E-16 0.0 ]
114 (%i8) transpose (u) . M - D . transpose (u);
115 [ 0.0 - 4.44089E-16 ]
120 @opencatbox{Categories:}
121 @category{Package lapack}
127 @deffn {Function} dgeqrf (@var{A})
129 Computes the QR decomposition of the matrix @var{A}.
130 All elements of @var{A} must be integer or floating point numbers.
131 @var{A} may or may not have the same number of rows and columns.
133 To make use of this function, you must load the LaPack package via
134 @code{load("lapack")}.
136 The real square matrix
137 m4_math(<<<\mathbf{A}>>>, <<<A>>>, <<<{\bf A}>>>)
141 <<<\mathbf{A} = \mathbf{Q}\mathbf{R}>>>,
143 <<<{\bf A} = {\bf Q} {\bf R}>>>)
146 m4_math(<<<{\bf Q}>>>, <<<Q>>>)
147 is a square orthogonal matrix with the same number of rows as
148 m4_math(<<<\mathbf{A}>>>, <<<A>>>, <<<{\bf A}>>>)
150 m4_math(<<<{\bf R}>>>, <<<R>>>)
151 is an upper triangular matrix and is the same size as
152 m4_mathdot(<<<{\bf A}>>>, A)
154 A list of two items is returned.
155 The first item is the matrix
156 m4_mathdot(<<<{\bf Q}>>>, Q)
157 The second item is the matrix
158 m4_mathcomma(<<<{\bf R}>>>, R)
159 The product @math{Q . R}, where "." is the noncommutative multiplication operator,
160 is equal to @var{A} (ignoring floating point round-off errors).
165 @c M : matrix ([1, -3.2, 8], [-11, 2.7, 5.9]);
166 @c [q, r] : dgeqrf (M);
171 (%i1) load ("lapack") $
172 (%i2) fpprintprec : 6 $
173 (%i3) M : matrix ([1, -3.2, 8], [-11, 2.7, 5.9]) $
174 (%i4) [q, r] : dgeqrf (M);
175 [ - .0905357 .995893 ]
178 [ - 11.0454 2.97863 5.15148 ]
180 [ 0 - 2.94241 8.50131 ]
182 [ - 7.77156E-16 1.77636E-15 - 8.88178E-16 ]
184 [ 0.0 - 1.33227E-15 8.88178E-16 ]
185 (%i6) mat_norm (%, 1);
189 @opencatbox{Categories:}
190 @category{Package lapack}
196 @deffn {Function} dgesv (@var{A}, @var{b})
198 Computes the solution @var{x} of the linear equation
199 m4_mathcomma(<<<{\bf A} x = b>>>, <<<@var{A} @var{x} = @var{b}>>>)
201 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
202 is a square matrix, and @math{b} is a matrix of the same number of rows
204 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
205 and any number of columns.
206 The return value @math{x} is the same size as @math{b}.
208 To make use of this function, you must load the LaPack package via
209 @code{load("lapack")}.
211 The elements of @var{A} and @var{b} must evaluate to real floating point numbers via @code{float};
212 thus elements may be any numeric type, symbolic numerical constants, or expressions which evaluate to floats.
213 The elements of @var{x} are always floating point numbers.
214 All arithmetic is carried out as floating point operations.
216 @code{dgesv} computes the solution via the LU decomposition of @var{A}.
220 @code{dgesv} computes the solution of the linear equation @math{@var{A} @var{x} = @var{b}}.
223 @c A : matrix ([1, -2.5], [0.375, 5]);
224 @c b : matrix ([1.75], [-0.625]);
226 @c dlange (inf_norm, b - A . x);
229 (%i1) A : matrix ([1, -2.5], [0.375, 5]);
233 (%i2) b : matrix ([1.75], [-0.625]);
237 (%i3) x : dgesv (A, b);
238 [ 1.210526315789474 ]
240 [ - 0.215789473684211 ]
241 (%i4) dlange (inf_norm, b - A.x);
245 @var{b} is a matrix with the same number of rows as @var{A} and any number of columns.
246 @var{x} is the same size as @var{b}.
249 @c A : matrix ([1, -0.15], [1.82, 2]);
250 @c b : matrix ([3.7, 1, 8], [-2.3, 5, -3.9]);
252 @c dlange (inf_norm, b - A . x);
255 (%i1) A : matrix ([1, -0.15], [1.82, 2]);
259 (%i2) b : matrix ([3.7, 1, 8], [-2.3, 5, -3.9]);
263 (%i3) x : dgesv (A, b);
264 [ 3.103827540695117 1.20985481742191 6.781786185657722 ]
266 [ -3.974483062032557 1.399032116146062 -8.121425428948527 ]
267 (%i4) dlange (inf_norm, b - A . x);
268 (%o4) 1.1102230246251565E-15
271 The elements of @var{A} and @var{b} must evaluate to real floating point numbers.
274 @c A : matrix ([5, -%pi], [1b0, 11/17]);
275 @c b : matrix ([%e], [sin(1)]);
277 @c dlange (inf_norm, b - A . x);
280 (%i1) A : matrix ([5, -%pi], [1b0, 11/17]);
286 (%i2) b : matrix ([%e], [sin(1)]);
290 (%i3) x : dgesv (A, b);
291 [ 0.690375643155986 ]
293 [ 0.233510982552952 ]
294 (%i4) dlange (inf_norm, b - A . x);
295 (%o4) 2.220446049250313E-16
298 @opencatbox{Categories:}
299 @category{Package lapack}
300 @category{Linear equations}
306 @deffn {Function} dgesvd (@var{A})
307 @deffnx {Function} dgesvd (@var{A}, @var{left_p}, @var{right_p})
309 Computes the singular value decomposition (SVD) of a matrix @var{A},
310 comprising the singular values and, optionally, the left and right singular vectors.
311 All elements of @var{A} must be integer or floating point numbers.
312 @var{A} might or might not be square (same number of rows and columns).
314 To make use of this function, you must load the LaPack package via
315 @code{load("lapack")}.
317 Let @math{m} be the number of rows, and @math{n} the number of columns of @var{A}.
318 The singular value decomposition of
319 m4_math(<<<\mathbf{A}>>>,<<<@math{A}>>>, <<<{\bf A}>>>)
320 comprises three matrices,
321 m4_mathcomma(<<<\mathbf{U}>>>,<<<U>>>,<<<{\bf U}>>>)
322 m4_mathcomma(<<<\mathbf{\Sigma}>>>,<<<@math{Sigma}>>>, <<<{\bf \Sigma}>>>)
324 m4_mathcomma(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
327 @c this code breaks texi2pdf: @math{@var{A} = @var{U} . @var{Sigma} . @var{V^T}}
328 @c @math{@var{A} = @var{U} . @var{Sigma} . @var{V}^T}
331 <<<\mathbf{A} = \mathbf{U} \mathbf{\Sigma} \mathbf{V}^T>>>,
332 <<<@math{{A} = {U} . {Sigma} . {V}^T}>>>,
333 <<<{\bf A} = {\bf U} {\bf \Sigma} {\bf V}^T>>>)
337 m4_math(<<<\mathbf{U}>>>, <<<@math{U}>>>, <<<{\bf U}>>>)
338 is an @math{m}-by-@math{m} unitary matrix,
339 m4_math(<<<\mathbf{\Sigma}>>>, <<<@math{Sigma}>>>, <<<{\bf\Sigma}>>>)
340 is an @math{m}-by-@math{n} diagonal matrix,
342 m4_math(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
343 is an @math{n}-by-@math{n} unitary matrix.
346 m4_math(<<<\mathbf{\sigma}_i>>>, <<<@math{sigma[i]}>>>, <<<{\bf \sigma}_i>>>)
347 be a diagonal element of
348 m4_mathcomma(<<<\mathbf{\Sigma}>>>, <<<@math{Sigma}>>>, <<<{\bf \Sigma}>>>)
350 m4_mathdot(<<<\mathbf{\Sigma}_{ii} = \sigma_i>>>,
351 <<<@math{@var{Sigma}[i, i] = @var{sigma}[i]}>>>,
352 <<<{\bf \Sigma}_{ii} = \sigma_i>>>)
354 m4_math(<<<\sigma_i>>>, <<<@math{sigma[i]}>>>)
355 are the so-called singular values of
356 m4_mathpunc(<<<;>>>, <<<\mathbf{A}>>>, <<<A>>>, <<<{\bf A}>>>)
357 these are real and nonnegative, and returned in descending order.
359 m4_math(<<<\min(m, n)>>>, <<<min(m, n)>>>)
361 m4_math(<<<\mathbf{U}>>>, <<<@math{U}>>>, <<<{\bf U}>>>)
363 m4_math(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
364 are the left and right singular vectors of
365 m4_mathdot(<<<\mathbf{A}>>>, <<<@math{A}>>>, <<<{\bf A}>>>)
366 Note that @code{dgesvd} returns the transpose of
367 m4_mathcomma(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
369 m4_math(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
372 @code{dgesvd(@var{A})} computes only the singular values of @var{A}.
373 @code{dgesvd(@var{A}, @var{left_p}, @var{right_p})} computes the singular values of @var{A}
374 and the left singular vectors when @math{@var{left_p} = @code{true}}
375 and the right singular vectors when @math{@var{right_p} = @code{true}}.
377 A list of three items is returned.
378 The first item is a list of the singular values.
379 The second item is @code{false} or the matrix of left singular vectors.
380 The third item is @code{false} or the matrix of right singular vectors.
387 @c M: matrix([1, 2, 3], [3.5, 0.5, 8], [-1, 2, -3], [4, 9, 7]);
389 @c [sigma, U, VT] : dgesvd (M, true, true);
393 @c genmatrix(lambda ([i, j], if i=j then sigma[i] else 0),
395 @c U . Sigma . VT - M;
396 @c transpose (U) . U;
397 @c VT . transpose (VT);
400 (%i1) load ("lapack")$
401 (%i2) fpprintprec : 6;
403 (%i3) M: matrix([1, 2, 3], [3.5, 0.5, 8], [-1, 2, -3], [4, 9, 7]);
412 (%o4) [[14.4744, 6.38637, .452547], false, false]
413 (%i5) [sigma, U, VT] : dgesvd (M, true, true);
414 (%o5) [[14.4744, 6.38637, .452547],
415 [ - .256731 .00816168 .959029 - .119523 ]
417 [ - .526456 .672116 - .206236 - .478091 ]
419 [ .107997 - .532278 - .0708315 - 0.83666 ]
421 [ - .803287 - .514659 - .180867 .239046 ]
422 [ - .374486 - .538209 - .755044 ]
424 [ .130623 - .836799 0.5317 ]]
426 [ - .917986 .100488 .383672 ]
427 (%i6) m : length (U);
429 (%i7) n : length (VT);
432 genmatrix(lambda ([i, j], if i=j then sigma[i] else 0),
441 (%i9) U . Sigma . VT - M;
442 [ 1.11022E-15 0.0 1.77636E-15 ]
444 [ 1.33227E-15 1.66533E-15 0.0 ]
446 [ - 4.44089E-16 - 8.88178E-16 4.44089E-16 ]
448 [ 8.88178E-16 1.77636E-15 8.88178E-16 ]
449 (%i10) transpose (U) . U;
450 [ 1.0 5.55112E-17 2.498E-16 2.77556E-17 ]
452 [ 5.55112E-17 1.0 5.55112E-17 4.16334E-17 ]
454 [ 2.498E-16 5.55112E-17 1.0 - 2.08167E-16 ]
456 [ 2.77556E-17 4.16334E-17 - 2.08167E-16 1.0 ]
457 (%i11) VT . transpose (VT);
458 [ 1.0 0.0 - 5.55112E-17 ]
460 (%o11) [ 0.0 1.0 5.55112E-17 ]
462 [ - 5.55112E-17 5.55112E-17 1.0 ]
465 @opencatbox{Categories:}
466 @category{Package lapack}
473 @deffn {Function} dlange (@var{norm}, @var{A})
474 @deffnx {Function} zlange (@var{norm}, @var{A})
476 Computes a norm or norm-like function of the matrix @var{A}. If
477 @var{A} is a real matrix, use @code{dlange}. For a matrix with
478 complex elements, use @code{zlange}.
480 To make use of this function, you must load the LaPack package via
481 @code{load("lapack")}.
483 @code{norm} specifies the kind of norm to be computed:
487 m4_math(<<<\max(|{\bf A}_{ij}|)>>>, <<<max(abs(A(i, j)))>>>)
488 where @math{i} and @math{j} range over
489 the rows and columns, respectively, of
490 m4_mathdot(<<<{\bf A}>>>, <<<A>>>)
491 Note that this function is not a proper matrix norm.
494 m4_math(<<<L_1>>>, <<<L[1]>>>)
496 m4_mathcomma(<<<{\bf A}>>>, <<<A>>>)
497 that is, the maximum of the sum of the absolute value of elements in each column.
500 m4_math(<<<L_\infty>>>, <<<L[inf]>>>)
502 m4_mathcomma(<<<{\bf A}>>>, <<<A>>>)
503 that is, the maximum of the sum of the absolute value of elements in each row.
505 Compute the Frobenius norm of
506 m4_mathcomma(<<<{\bf A}>>>, <<<A>>>)
507 that is, the square root of the sum of squares of the matrix elements.
510 @opencatbox{Categories:}
511 @category{Package lapack}
517 @deffn {Function} dgemm (@var{A}, @var{B})
518 @deffnx {Function} dgemm (@var{A}, @var{B}, @var{options})
519 Compute the product of two matrices and optionally add the product to
522 In the simplest form, @code{dgemm(@var{A}, @var{B})} computes the
523 product of the two real matrices, @var{A} and @var{B}.
525 To make use of this function, you must load the LaPack package via
526 @code{load("lapack")}.
528 In the second form, @code{dgemm} computes
529 m4_math(<<<\alpha {\bf A} {\bf B} + \beta {\bf C}>>>,
530 <<<@var{alpha} * @var{A} * @var{B} + @var{beta} * @var{C}>>>)
532 m4_mathcomma(<<<{\bf A}>>>, <<<@var{A}>>>)
533 m4_mathcomma(<<<{\bf B}>>>, <<<@var{B}>>>)
535 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
536 are real matrices of the appropriate sizes and
537 m4_math(<<<\alpha>>>, <<<alpha>>>)
539 m4_math(<<<\beta>>>, <<<beta>>>)
540 are real numbers. Optionally,
541 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
543 m4_math(<<<{\bf B}>>>, <<<@var{B}>>>)
545 be transposed before computing the product. The extra parameters are
546 specified by optional keyword arguments: The keyword arguments are
547 optional and may be specified in any order. They all take the form
548 @code{key=val}. The keyword arguments are:
553 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
554 that should be added. The default is @code{false},
555 which means no matrix is added.
558 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
560 m4_math(<<<{\bf B}>>>, <<<@var{B}>>>)
561 is multiplied by this value. The
565 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
566 is given, this value multiplies
567 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
569 is added. The default value is 0, which implies that
570 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
573 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
574 is given. Hence, be sure to specify a non-zero
576 m4_mathdot(<<<\beta>>>, <<<@math{beta}>>>)
578 If @code{true}, the transpose of
579 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
581 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
582 for the product. The default is @code{false}.
584 If @code{true}, the transpose of
585 m4_math(<<<{\bf B}>>>, <<<@var{B}>>>)
587 m4_math(<<<{\bf B}>>>, <<<@var{B}>>>)
588 for the product. The default is @code{false}.
593 @c A : matrix([1,2,3],[4,5,6],[7,8,9]);
594 @c B : matrix([-1,-2,-3],[-4,-5,-6],[-7,-8,-9]);
595 @c C : matrix([3,2,1],[6,5,4],[9,8,7]);
598 @c dgemm(A,B,transpose_a=true);
600 @c dgemm(A,B,c=C,beta=1);
602 @c dgemm(A,B,c=C,beta=1, alpha=-1);
606 (%i1) load ("lapack")$
607 (%i2) A : matrix([1,2,3],[4,5,6],[7,8,9]);
613 (%i3) B : matrix([-1,-2,-3],[-4,-5,-6],[-7,-8,-9]);
616 (%o3) [ - 4 - 5 - 6 ]
619 (%i4) C : matrix([3,2,1],[6,5,4],[9,8,7]);
626 [ - 30.0 - 36.0 - 42.0 ]
628 (%o5) [ - 66.0 - 81.0 - 96.0 ]
630 [ - 102.0 - 126.0 - 150.0 ]
634 (%o6) [ - 66 - 81 - 96 ]
636 [ - 102 - 126 - 150 ]
637 (%i7) dgemm(A,B,transpose_a=true);
638 [ - 66.0 - 78.0 - 90.0 ]
640 (%o7) [ - 78.0 - 93.0 - 108.0 ]
642 [ - 90.0 - 108.0 - 126.0 ]
643 (%i8) transpose(A) . B;
646 (%o8) [ - 78 - 93 - 108 ]
649 (%i9) dgemm(A,B,c=C,beta=1);
650 [ - 27.0 - 34.0 - 41.0 ]
652 (%o9) [ - 60.0 - 76.0 - 92.0 ]
654 [ - 93.0 - 118.0 - 143.0 ]
658 (%o10) [ - 60 - 76 - 92 ]
661 (%i11) dgemm(A,B,c=C,beta=1, alpha=-1);
664 (%o11) [ 72.0 86.0 100.0 ]
666 [ 111.0 134.0 157.0 ]
675 @opencatbox{Categories:}
676 @category{Package lapack}
682 @deffn {Function} zgeev (@var{A})
683 @deffnx {Function} zgeev (@var{A}, @var{right_p}, @var{left_p})
685 Like @mrefcomma{dgeev} but the matrix
686 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
689 To make use of this function, you must load the LaPack package via
690 @code{load("lapack")}.
692 @opencatbox{Categories:}
693 @category{Package lapack}
699 @deffn {Function} zheev (@var{A})
700 @deffnx {Function} zheev (@var{A}, @var{eigvec_p})
702 Like @mrefcomma{dgeev} but the matrix
703 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
704 is assumed to be a square
705 complex Hermitian matrix. If @var{eigvec_p} is @code{true}, then the
706 eigenvectors of the matrix are also computed.
708 To make use of this function, you must load the LaPack package via
709 @code{load("lapack")}.
711 No check is made that the matrix
712 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
713 is, in fact, Hermitian.
715 A list of two items is returned, as in @code{dgeev}: a list of
716 eigenvalues, and @code{false} or the matrix of the eigenvectors.
721 @c [9.14 +%i*0.00 , -4.37 -%i*9.22 , -1.98 -%i*1.72 , -8.96 -%i*9.50],
722 @c [-4.37 +%i*9.22 , -3.35 +%i*0.00 , 2.25 -%i*9.51 , 2.57 +%i*2.40],
723 @c [-1.98 +%i*1.72 , 2.25 +%i*9.51 , -4.82 +%i*0.00 , -3.24 +%i*2.04],
724 @c [-8.96 +%i*9.50 , 2.57 -%i*2.40 , -3.24 -%i*2.04 , 8.44 +%i*0.00]);
730 An example of computing the eigenvalues and then eigenvalues and
731 eigenvectors of an Hermitian matrix.
733 (%i1) load("lapack")$
735 [9.14 +%i*0.00 , -4.37 -%i*9.22 , -1.98 -%i*1.72 , -8.96 -%i*9.50],
736 [-4.37 +%i*9.22 , -3.35 +%i*0.00 , 2.25 -%i*9.51 , 2.57 +%i*2.40],
737 [-1.98 +%i*1.72 , 2.25 +%i*9.51 , -4.82 +%i*0.00 , -3.24 +%i*2.04],
738 [-8.96 +%i*9.50 , 2.57 -%i*2.40 , -3.24 -%i*2.04 , 8.44 +%i*0.00]);
740 [ 9.14 (- 9.22 %i) - 4.37 (- 1.72 %i) - 1.98 (- 9.5 %i) - 8.96 ]
742 [ 9.22 %i - 4.37 - 3.35 2.25 - 9.51 %i 2.4 %i + 2.57 ]
744 [ 1.72 %i - 1.98 9.51 %i + 2.25 - 4.82 2.04 %i - 3.24 ]
746 [ 9.5 %i - 8.96 2.57 - 2.4 %i (- 2.04 %i) - 3.24 8.44 ]
748 (%o3) [[- 16.00474647209473, - 6.764970154793324, 6.665711453507098,
749 25.51400517338097], false]
750 (%i4) E:zheev(A,true)$
752 (%o5) [- 16.00474647209474, - 6.764970154793325, 6.665711453507101,
755 [ 0.2674650533172745 %i + 0.2175453586665017 ]
757 [ 0.002696730886619885 %i + 0.6968836773391712 ]
759 [ (- 0.6082406376714117 %i) - 0.01210614292697931 ]
761 [ 0.1593081858095037 ]
762 [ 0.2644937470667444 %i + 0.4773693349937472 ]
764 [ (- 0.2852389036031621 %i) - 0.1414362742011673 ]
766 [ 0.2654607680986639 %i + 0.4467818117184174 ]
768 [ 0.5750762708542709 ]
769 [ 0.2810649767305922 %i - 0.1335263928245182 ]
771 [ 0.2866310132869556 %i - 0.4536971347853274 ]
773 [ (- 0.2933684323754295 %i) - 0.4954972425541057 ]
775 [ 0.5325337537576771 ]
776 [ (- 0.5737316575503476 %i) - 0.3966146799427706 ]
778 [ 0.01826502619021457 %i + 0.3530557704387017 ]
780 [ 0.1673700900085425 %i + 0.01476684746229564 ]
782 [ 0.6002632636961784 ]
785 @opencatbox{Categories:}
786 @category{Package lapack}