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 @
24 @fname{dgeev} (@var{A}) @
25 @fname{dgeev} (@var{A}, @var{right_p}, @var{left_p})
27 Computes the eigenvalues and, optionally, the eigenvectors of a matrix @var{A}.
28 All elements of @var{A} must be integer or floating point numbers.
29 @var{A} must be square (same number of rows and columns).
30 @var{A} might or might not be symmetric.
32 To make use of this function, you must load the LaPack package via
33 @code{load("lapack")}.
35 @code{dgeev(@var{A})} computes only the eigenvalues of @var{A}.
36 @code{dgeev(@var{A}, @var{right_p}, @var{left_p})} computes the eigenvalues of @var{A}
37 and the right eigenvectors when @math{@var{right_p} = @code{true}}
38 and the left eigenvectors when @math{@var{left_p} = @code{true}}.
40 A list of three items is returned.
41 The first item is a list of the eigenvalues.
42 The second item is @code{false} or the matrix of right eigenvectors.
43 The third item is @code{false} or the matrix of left eigenvectors.
46 m4_math(<<<v_j>>>,<<<v_j>>>)
47 (the @math{j}-th column of the right eigenvector matrix) satisfies
50 <<<\mathbf{A} v_j = \lambda_j v_j>>>,
51 <<<@math{A . v_j = lambda_j . v_j}>>>,
52 <<<{\bf A} v_j = \lambda_j v_j>>>
56 m4_math(<<<\lambda_j>>>,<<<@math{lambda_j}>>>)
57 is the corresponding eigenvalue.
59 m4_math(<<<u_j>>>,<<<u_j>>>)
60 (the @math{j}-th column of the left eigenvector matrix) satisfies
63 <<<u_j^\mathbf{H} \mathbf{A} = \lambda_j u_j^\mathbf{H}>>>,
64 <<<@math{u_j^H . A = lambda_j . u_j^H}>>>,
65 <<<u_j^{\bf H} {\bf A} = \lambda_j u_j^{\bf H}>>>
69 m4_math(<<<u_j^\mathbf{H}>>>,<<<u_j^H>>>, <<<u_j^{\bf H}>>>)
70 denotes the conjugate transpose of
71 m4_mathdot(<<<u_j>>>,<<<u_j>>>)
72 The Maxima function @code{ctranspose} computes the conjugate transpose.
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}>>>)
140 <<<\mathbf{A} = \mathbf{Q}\mathbf{R}>>>,
142 <<<{\bf A} = {\bf Q} {\bf R}>>>)
145 m4_math(<<<{\bf Q}>>>, {{{Q}}})
146 is a square orthogonal matrix with the same number of rows as
147 m4_math(<<<\mathbf{A}>>>, <<<A>>>, <<<{\bf A}>>>)
149 m4_math(<<<{\bf R}>>>, {{{R}}})
150 is an upper triangular matrix.
152 A list of two items is returned.
153 The first item is the matrix @var{Q}.
154 The second item is the matrix @var{R}, which is the same size as @var{A},
155 and which has all elements equal to zero below the diagonal.
156 The product @code{@var{Q} . @var{R}}, where "." is the noncommutative multiplication operator,
157 is equal to @var{A} (ignoring floating point round-off errors).
162 @c M : matrix ([1, -3.2, 8], [-11, 2.7, 5.9]);
163 @c [q, r] : dgeqrf (M);
168 (%i1) load ("lapack") $
169 (%i2) fpprintprec : 6 $
170 (%i3) M : matrix ([1, -3.2, 8], [-11, 2.7, 5.9]) $
171 (%i4) [q, r] : dgeqrf (M);
172 [ - .0905357 .995893 ]
175 [ - 11.0454 2.97863 5.15148 ]
177 [ 0 - 2.94241 8.50131 ]
179 [ - 7.77156E-16 1.77636E-15 - 8.88178E-16 ]
181 [ 0.0 - 1.33227E-15 8.88178E-16 ]
182 (%i6) mat_norm (%, 1);
186 @opencatbox{Categories:}
187 @category{Package lapack}
193 @deffn {Function} dgesv (@var{A}, @var{b})
195 Computes the solution @var{x} of the linear equation @math{@var{A} @var{x} = @var{b}},
196 where @var{A} is a square matrix, and @var{b} is a matrix of the same number of rows
197 as @var{A} and any number of columns.
198 The return value @var{x} is the same size as @var{b}.
200 To make use of this function, you must load the LaPack package via
201 @code{load("lapack")}.
203 The elements of @var{A} and @var{b} must evaluate to real floating point numbers via @code{float};
204 thus elements may be any numeric type, symbolic numerical constants, or expressions which evaluate to floats.
205 The elements of @var{x} are always floating point numbers.
206 All arithmetic is carried out as floating point operations.
208 @code{dgesv} computes the solution via the LU decomposition of @var{A}.
212 @code{dgesv} computes the solution of the linear equation @math{@var{A} @var{x} = @var{b}}.
215 @c A : matrix ([1, -2.5], [0.375, 5]);
216 @c b : matrix ([1.75], [-0.625]);
218 @c dlange (inf_norm, b - A . x);
221 (%i1) A : matrix ([1, -2.5], [0.375, 5]);
225 (%i2) b : matrix ([1.75], [-0.625]);
229 (%i3) x : dgesv (A, b);
230 [ 1.210526315789474 ]
232 [ - 0.215789473684211 ]
233 (%i4) dlange (inf_norm, b - A.x);
237 @var{b} is a matrix with the same number of rows as @var{A} and any number of columns.
238 @var{x} is the same size as @var{b}.
241 @c A : matrix ([1, -0.15], [1.82, 2]);
242 @c b : matrix ([3.7, 1, 8], [-2.3, 5, -3.9]);
244 @c dlange (inf_norm, b - A . x);
247 (%i1) A : matrix ([1, -0.15], [1.82, 2]);
251 (%i2) b : matrix ([3.7, 1, 8], [-2.3, 5, -3.9]);
255 (%i3) x : dgesv (A, b);
256 [ 3.103827540695117 1.20985481742191 6.781786185657722 ]
258 [ -3.974483062032557 1.399032116146062 -8.121425428948527 ]
259 (%i4) dlange (inf_norm, b - A . x);
260 (%o4) 1.1102230246251565E-15
263 The elements of @var{A} and @var{b} must evaluate to real floating point numbers.
266 @c A : matrix ([5, -%pi], [1b0, 11/17]);
267 @c b : matrix ([%e], [sin(1)]);
269 @c dlange (inf_norm, b - A . x);
272 (%i1) A : matrix ([5, -%pi], [1b0, 11/17]);
278 (%i2) b : matrix ([%e], [sin(1)]);
282 (%i3) x : dgesv (A, b);
283 [ 0.690375643155986 ]
285 [ 0.233510982552952 ]
286 (%i4) dlange (inf_norm, b - A . x);
287 (%o4) 2.220446049250313E-16
290 @opencatbox{Categories:}
291 @category{Package lapack}
292 @category{Linear equations}
298 @deffn {Function} dgesvd @
299 @fname{dgesvd} (@var{A}) @
300 @fname{dgesvd} (@var{A}, @var{left_p}, @var{right_p})
302 Computes the singular value decomposition (SVD) of a matrix @var{A},
303 comprising the singular values and, optionally, the left and right singular vectors.
304 All elements of @var{A} must be integer or floating point numbers.
305 @var{A} might or might not be square (same number of rows and columns).
307 To make use of this function, you must load the LaPack package via
308 @code{load("lapack")}.
310 Let @math{m} be the number of rows, and @math{n} the number of columns of @var{A}.
311 The singular value decomposition of
312 m4_math(<<<\mathbf{A}>>>,<<<@math{A}>>>, <<<{\bf A}>>>)
313 comprises three matrices,
314 m4_mathcomma(<<<\mathbf{U}>>>,<<<U>>>,<<<{\bf U}>>>)
315 m4_mathcomma(<<<\mathbf{\Sigma}>>>,<<<@math{Sigma}>>>, <<<{\bf \Sigma}>>>)
317 m4_mathcomma(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
320 @c this code breaks texi2pdf: @math{@var{A} = @var{U} . @var{Sigma} . @var{V^T}}
321 @c @math{@var{A} = @var{U} . @var{Sigma} . @var{V}^T}
324 <<<\mathbf{A} = \mathbf{U} \mathbf{\Sigma} \mathbf{V}^T>>>,
325 <<<@math{@var{A} = @var{U} . @var{Sigma} . @var{V}^T}>>>,
326 <<<{\bf A} = {\bf U} {\bf \Sigma} {\bf V}^T>>>)
330 m4_math(<<<\mathbf{U}>>>, <<<@math{U}>>>, <<<{\bf U}>>>)
331 is an @math{m}-by-@math{m} unitary matrix,
332 m4_math(<<<\mathbf{\Sigma}>>>, <<<@math{Sigma}>>>, <<<{\bf\Sigma}>>>)
333 is an @math{m}-by-@math{n} diagonal matrix,
335 m4_math(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
336 is an @math{n}-by-@math{n} unitary matrix.
339 m4_math(<<<\mathbf{\sigma}_i>>>, <<<@math{sigma[i]}>>>, <<<{\bf \sigma}_i>>>)
340 be a diagonal element of
341 m4_mathcomma(<<<\mathbf{\Sigma}>>>, <<<@math{Sigma}>>>, <<<{\bf \Sigma}>>>)
343 m4_mathdot(<<<\mathbf{\Sigma}_{ii} = \sigma_i>>>,
344 <<<@math{@var{Sigma}[i, i] = @var{sigma}[i]}>>>,
345 <<<{\bf \Sigma}_{ii} = \sigma_i>>>)
347 m4_math(<<<\sigma_i>>>, <<<@math{sigma[i]}>>>)
348 are the so-called singular values of
349 m4_mathpunc(<<<;>>>, <<<\mathbf{A}>>>, <<<A>>>, <<<{\bf A}>>>)
350 these are real and nonnegative, and returned in descending order.
352 m4_math(<<<\min(m, n)>>>, <<<min(m, n)>>>)
354 m4_math(<<<\mathbf{U}>>>, <<<@math{U}>>>, <<<{\bf U}>>>)
356 m4_math(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
357 are the left and right singular vectors of
358 m4_mathdot(<<<\mathbf{A}>>>, <<<@math{A}>>>, <<<{\bf A}>>>)
359 Note that @code{dgesvd} returns the transpose of
360 m4_mathcomma(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
362 m4_math(<<<\mathbf{V}>>>, <<<@math{V}>>>, <<<{\bf V}>>>)
365 @code{dgesvd(@var{A})} computes only the singular values of @var{A}.
366 @code{dgesvd(@var{A}, @var{left_p}, @var{right_p})} computes the singular values of @var{A}
367 and the left singular vectors when @math{@var{left_p} = @code{true}}
368 and the right singular vectors when @math{@var{right_p} = @code{true}}.
370 A list of three items is returned.
371 The first item is a list of the singular values.
372 The second item is @code{false} or the matrix of left singular vectors.
373 The third item is @code{false} or the matrix of right singular vectors.
380 @c M: matrix([1, 2, 3], [3.5, 0.5, 8], [-1, 2, -3], [4, 9, 7]);
382 @c [sigma, U, VT] : dgesvd (M, true, true);
386 @c genmatrix(lambda ([i, j], if i=j then sigma[i] else 0),
388 @c U . Sigma . VT - M;
389 @c transpose (U) . U;
390 @c VT . transpose (VT);
393 (%i1) load ("lapack")$
394 (%i2) fpprintprec : 6;
396 (%i3) M: matrix([1, 2, 3], [3.5, 0.5, 8], [-1, 2, -3], [4, 9, 7]);
405 (%o4) [[14.4744, 6.38637, .452547], false, false]
406 (%i5) [sigma, U, VT] : dgesvd (M, true, true);
407 (%o5) [[14.4744, 6.38637, .452547],
408 [ - .256731 .00816168 .959029 - .119523 ]
410 [ - .526456 .672116 - .206236 - .478091 ]
412 [ .107997 - .532278 - .0708315 - 0.83666 ]
414 [ - .803287 - .514659 - .180867 .239046 ]
415 [ - .374486 - .538209 - .755044 ]
417 [ .130623 - .836799 0.5317 ]]
419 [ - .917986 .100488 .383672 ]
420 (%i6) m : length (U);
422 (%i7) n : length (VT);
425 genmatrix(lambda ([i, j], if i=j then sigma[i] else 0),
434 (%i9) U . Sigma . VT - M;
435 [ 1.11022E-15 0.0 1.77636E-15 ]
437 [ 1.33227E-15 1.66533E-15 0.0 ]
439 [ - 4.44089E-16 - 8.88178E-16 4.44089E-16 ]
441 [ 8.88178E-16 1.77636E-15 8.88178E-16 ]
442 (%i10) transpose (U) . U;
443 [ 1.0 5.55112E-17 2.498E-16 2.77556E-17 ]
445 [ 5.55112E-17 1.0 5.55112E-17 4.16334E-17 ]
447 [ 2.498E-16 5.55112E-17 1.0 - 2.08167E-16 ]
449 [ 2.77556E-17 4.16334E-17 - 2.08167E-16 1.0 ]
450 (%i11) VT . transpose (VT);
451 [ 1.0 0.0 - 5.55112E-17 ]
453 (%o11) [ 0.0 1.0 5.55112E-17 ]
455 [ - 5.55112E-17 5.55112E-17 1.0 ]
458 @opencatbox{Categories:}
459 @category{Package lapack}
466 @deffn {Function} dlange (@var{norm}, @var{A})
467 @deffnx {Function} zlange (@var{norm}, @var{A})
469 Computes a norm or norm-like function of the matrix @var{A}. If
470 @var{A} is a real matrix, use @code{dlange}. For a matrix with
471 complex elements, use @code{zlange}.
473 To make use of this function, you must load the LaPack package via
474 @code{load("lapack")}.
476 @code{norm} specifies the kind of norm to be computed:
479 Compute @math{max(abs(A(i, j)))} where @math{i} and @math{j} range over
480 the rows and columns, respectively, of @var{A}.
481 Note that this function is not a proper matrix norm.
483 Compute the @math{L[1]} norm of @var{A},
484 that is, the maximum of the sum of the absolute value of elements in each column.
486 Compute the @math{L[inf]} norm of @var{A},
487 that is, the maximum of the sum of the absolute value of elements in each row.
489 Compute the Frobenius norm of @var{A},
490 that is, the square root of the sum of squares of the matrix elements.
493 @opencatbox{Categories:}
494 @category{Package lapack}
500 @deffn {Function} dgemm @
501 @fname{dgemm} (@var{A}, @var{B}) @
502 @fname{dgemm} (@var{A}, @var{B}, @var{options})
503 Compute the product of two matrices and optionally add the product to
506 In the simplest form, @code{dgemm(@var{A}, @var{B})} computes the
507 product of the two real matrices, @var{A} and @var{B}.
509 To make use of this function, you must load the LaPack package via
510 @code{load("lapack")}.
512 In the second form, @code{dgemm} computes the @math{@var{alpha} *
513 @var{A} * @var{B} + @var{beta} * @var{C}} where @var{A}, @var{B},
514 @var{C} are real matrices of the appropriate sizes and @var{alpha} and
515 @var{beta} are real numbers. Optionally, @var{A} and/or @var{B} can
516 be transposed before computing the product. The extra parameters are
517 specified by optional keyword arguments: The keyword arguments are
518 optional and may be specified in any order. They all take the form
519 @code{key=val}. The keyword arguments are:
523 The matrix @var{C} that should be added. The default is @code{false},
524 which means no matrix is added.
526 The product of @var{A} and @var{B} is multiplied by this value. The
529 If a matrix @var{C} is given, this value multiplies @var{C} before it
530 is added. The default value is 0, which implies that @var{C} is not
531 added, even if @var{C} is given. Hence, be sure to specify a non-zero
532 value for @var{beta}.
534 If @code{true}, the transpose of @var{A} is used instead of @var{A}
535 for the product. The default is @code{false}.
537 If @code{true}, the transpose of @var{B} is used instead of @var{B}
538 for the product. The default is @code{false}.
543 @c A : matrix([1,2,3],[4,5,6],[7,8,9]);
544 @c B : matrix([-1,-2,-3],[-4,-5,-6],[-7,-8,-9]);
545 @c C : matrix([3,2,1],[6,5,4],[9,8,7]);
548 @c dgemm(A,B,transpose_a=true);
550 @c dgemm(A,B,c=C,beta=1);
552 @c dgemm(A,B,c=C,beta=1, alpha=-1);
556 (%i1) load ("lapack")$
557 (%i2) A : matrix([1,2,3],[4,5,6],[7,8,9]);
563 (%i3) B : matrix([-1,-2,-3],[-4,-5,-6],[-7,-8,-9]);
566 (%o3) [ - 4 - 5 - 6 ]
569 (%i4) C : matrix([3,2,1],[6,5,4],[9,8,7]);
576 [ - 30.0 - 36.0 - 42.0 ]
578 (%o5) [ - 66.0 - 81.0 - 96.0 ]
580 [ - 102.0 - 126.0 - 150.0 ]
584 (%o6) [ - 66 - 81 - 96 ]
586 [ - 102 - 126 - 150 ]
587 (%i7) dgemm(A,B,transpose_a=true);
588 [ - 66.0 - 78.0 - 90.0 ]
590 (%o7) [ - 78.0 - 93.0 - 108.0 ]
592 [ - 90.0 - 108.0 - 126.0 ]
593 (%i8) transpose(A) . B;
596 (%o8) [ - 78 - 93 - 108 ]
599 (%i9) dgemm(A,B,c=C,beta=1);
600 [ - 27.0 - 34.0 - 41.0 ]
602 (%o9) [ - 60.0 - 76.0 - 92.0 ]
604 [ - 93.0 - 118.0 - 143.0 ]
608 (%o10) [ - 60 - 76 - 92 ]
611 (%i11) dgemm(A,B,c=C,beta=1, alpha=-1);
614 (%o11) [ 72.0 86.0 100.0 ]
616 [ 111.0 134.0 157.0 ]
625 @opencatbox{Categories:}
626 @category{Package lapack}
632 @deffn {Function} zgeev @
633 @fname{zgeev} (@var{A}) @
634 @fname{zgeev} (@var{A}, @var{right_p}, @var{left_p})
636 Like @code{dgeev}, but the matrix @var{A} is complex.
638 To make use of this function, you must load the LaPack package via
639 @code{load("lapack")}.
641 @opencatbox{Categories:}
642 @category{Package lapack}
648 @deffn {Function} zheev @
649 @fname{zheev} (@var{A}) @
650 @fname{zheev} (@var{A}, @var{eigvec_p})
652 Like @code{dgeev}, but the matrix @var{A} is assumed to be a square
653 complex Hermitian matrix. If @var{eigvec_p} is @code{true}, then the
654 eigenvectors of the matrix are also computed.
656 To make use of this function, you must load the LaPack package via
657 @code{load("lapack")}.
659 No check is made that the matrix @var{A} is, in fact, Hermitian.
661 A list of two items is returned, as in @code{dgeev}: a list of
662 eigenvalues, and @code{false} or the matrix of the eigenvectors.
667 @c [9.14 +%i*0.00 , -4.37 -%i*9.22 , -1.98 -%i*1.72 , -8.96 -%i*9.50],
668 @c [-4.37 +%i*9.22 , -3.35 +%i*0.00 , 2.25 -%i*9.51 , 2.57 +%i*2.40],
669 @c [-1.98 +%i*1.72 , 2.25 +%i*9.51 , -4.82 +%i*0.00 , -3.24 +%i*2.04],
670 @c [-8.96 +%i*9.50 , 2.57 -%i*2.40 , -3.24 -%i*2.04 , 8.44 +%i*0.00]);
675 An example of computing the eigenvalues and then eigenvalues and
676 eigenvectors of an Hermitian matrix.
678 (%i1) load("lapack")$
680 [9.14 +%i*0.00 , -4.37 -%i*9.22 , -1.98 -%i*1.72 , -8.96 -%i*9.50],
681 [-4.37 +%i*9.22 , -3.35 +%i*0.00 , 2.25 -%i*9.51 , 2.57 +%i*2.40],
682 [-1.98 +%i*1.72 , 2.25 +%i*9.51 , -4.82 +%i*0.00 , -3.24 +%i*2.04],
683 [-8.96 +%i*9.50 , 2.57 -%i*2.40 , -3.24 -%i*2.04 , 8.44 +%i*0.00]);
685 [ 9.14 (- 9.22 %i) - 4.37 (- 1.72 %i) - 1.98 (- 9.5 %i) - 8.96 ]
687 [ 9.22 %i - 4.37 - 3.35 2.25 - 9.51 %i 2.4 %i + 2.57 ]
689 [ 1.72 %i - 1.98 9.51 %i + 2.25 - 4.82 2.04 %i - 3.24 ]
691 [ 9.5 %i - 8.96 2.57 - 2.4 %i (- 2.04 %i) - 3.24 8.44 ]
693 (%o3) [[- 16.00474647209473, - 6.764970154793324, 6.665711453507098,
694 25.51400517338097], false]
695 (%i4) E:zheev(A,true)$
697 (%o5) [- 16.00474647209474, - 6.764970154793325, 6.665711453507101,
700 [ 0.2674650533172745 %i + 0.2175453586665017 ]
702 [ 0.002696730886619885 %i + 0.6968836773391712 ]
704 [ (- 0.6082406376714117 %i) - 0.01210614292697931 ]
706 [ 0.1593081858095037 ]
707 [ 0.2644937470667444 %i + 0.4773693349937472 ]
709 [ (- 0.2852389036031621 %i) - 0.1414362742011673 ]
711 [ 0.2654607680986639 %i + 0.4467818117184174 ]
713 [ 0.5750762708542709 ]
714 [ 0.2810649767305922 %i - 0.1335263928245182 ]
716 [ 0.2866310132869556 %i - 0.4536971347853274 ]
718 [ (- 0.2933684323754295 %i) - 0.4954972425541057 ]
720 [ 0.5325337537576771 ]
721 [ (- 0.5737316575503476 %i) - 0.3966146799427706 ]
723 [ 0.01826502619021457 %i + 0.3530557704387017 ]
725 [ 0.1673700900085425 %i + 0.01476684746229564 ]
727 [ 0.6002632636961784 ]
730 @opencatbox{Categories:}
731 @category{Package lapack}