Revise and extend documentation for function trace_options.
[maxima.git] / doc / info / lapack.texi.m4
blob88f45d9b8c64c593a594d063ec78eaa5ba913449
1 @menu
2 * Introduction to lapack::
3 * Functions and Variables for lapack::
4 @end menu
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}
16 @closecatbox
19 @node Functions and Variables for lapack, , Introduction to lapack, Package lapack
20 @section Functions and Variables for lapack
22 @anchor{dgeev}
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.
44 The right eigenvector
45 m4_math(<<<v_j>>>,<<<v_j>>>)
46 (the @math{j}-th column of the right eigenvector matrix) satisfies
48 m4_displaymath(
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>>>
54 where
55 m4_math(<<<\lambda_j>>>,<<<@math{lambda_j}>>>)
56 is the corresponding eigenvalue.
57 The left eigenvector
58 m4_math(<<<u_j>>>,<<<u_j>>>)
59 (the @math{j}-th column of the left eigenvector matrix) satisfies
61 m4_displaymath(
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}>>>
67 where
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.
77 Example:
79 @c ===beg===
80 @c load ("lapack")$
81 @c fpprintprec : 6;
82 @c M : matrix ([9.5, 1.75], [3.25, 10.45]);
83 @c dgeev (M);
84 @c [L, v, u] : dgeev (M, true, true);
85 @c D : apply (diag_matrix, L);
86 @c M . v - v . D;
87 @c transpose (u) . M - D . transpose (u);
88 @c ===end===
89 @example
90 (%i1) load ("lapack")$
91 @group
92 (%i2) fpprintprec : 6;
93 (%o2)                           6
94 @end group
95 @group
96 (%i3) M : matrix ([9.5, 1.75], [3.25, 10.45]);
97                          [ 9.5   1.75  ]
98 (%o3)                    [             ]
99                          [ 3.25  10.45 ]
100 @end group
101 @group
102 (%i4) dgeev (M);
103 (%o4)          [[7.54331, 12.4067], false, false]
104 @end group
105 @group
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 ]
111                                       [                        ]]
112                                       [  0.515792   - 0.666642 ]
113 @end group
114 @group
115 (%i6) D : apply (diag_matrix, L);
116                       [ 7.54331     0    ]
117 (%o6)                 [                  ]
118                       [    0     12.4067 ]
119 @end group
120 @group
121 (%i7) M . v - v . D;
122                 [      0.0       - 8.88178e-16 ]
123 (%o7)           [                              ]
124                 [ - 8.88178e-16       0.0      ]
125 @end group
126 @group
127 (%i8) transpose (u) . M - D . transpose (u);
128                      [ 0.0  - 4.44089e-16 ]
129 (%o8)                [                    ]
130                      [ 0.0       0.0      ]
131 @end group
132 @end example
134 @opencatbox{Categories:}
135 @category{Package lapack}
136 @closecatbox
138 @end deffn
140 @anchor{dgeqrf}
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}>>>)
152 can be decomposed as
154 m4_displaymath(
155 <<<\mathbf{A} = \mathbf{Q}\mathbf{R}>>>,
156 <<<A = QR>>>,
157 <<<{\bf A} = {\bf Q} {\bf R}>>>)
159 where
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).
176 @c ===beg===
177 @c load ("lapack")$
178 @c fpprintprec : 6;
179 @c M : matrix ([1, -3.2, 8], [-11, 2.7, 5.9]);
180 @c [q, r] : dgeqrf (M);
181 @c q . r - M;
182 @c mat_norm (%, 1);
183 @c ===end===
184 @example
185 (%i1) load ("lapack")$
186 @group
187 (%i2) fpprintprec : 6;
188 (%o2)                           6
189 @end group
190 @group
191 (%i3) M : matrix ([1, -3.2, 8], [-11, 2.7, 5.9]);
192                       [  1    - 3.2   8  ]
193 (%o3)                 [                  ]
194                       [ - 11   2.7   5.9 ]
195 @end group
196 @group
197 (%i4) [q, r] : dgeqrf (M);
198        [ - 0.0905357  0.995893  ]
199 (%o4) [[                        ], 
200        [  0.995893    0.0905357 ]
201                                [ - 11.0454   2.97863   5.15148 ]
202                                [                               ]]
203                                [    0.0     - 2.94241  8.50131 ]
204 @end group
205 @group
206 (%i5) q . r - M;
207          [ - 7.77156e-16   1.77636e-15   - 8.88178e-16 ]
208 (%o5)    [                                             ]
209          [      0.0       - 1.33227e-15   8.88178e-16  ]
210 @end group
211 @group
212 (%i6) mat_norm (%, 1);
213 (%o6)                      3.10862e-15
214 @end group
215 @end example
217 @opencatbox{Categories:}
218 @category{Package lapack}
219 @closecatbox
221 @end deffn
223 @anchor{dgesv}
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}>>>)
228 where
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}.
246 Examples:
248 @code{dgesv} computes the solution of the linear equation @math{@var{A} @var{x} = @var{b}}.
250 @c ===beg===
251 @c A : matrix ([1, -2.5], [0.375, 5]);
252 @c b : matrix ([1.75], [-0.625]);
253 @c x : dgesv (A, b);
254 @c dlange (inf_norm, b - A . x);
255 @c ===end===
256 @example
257 @group
258 (%i1) A : matrix ([1, -2.5], [0.375, 5]);
259                         [   1    - 2.5 ]
260 (%o1)                   [              ]
261                         [ 0.375    5   ]
262 @end group
263 @group
264 (%i2) b : matrix ([1.75], [-0.625]);
265                            [  1.75   ]
266 (%o2)                      [         ]
267                            [ - 0.625 ]
268 @end group
269 @group
270 (%i3) x : dgesv (A, b);
271                     [   1    - 2.5 ]  [  1.75   ]
272 (%o3)         dgesv([              ], [         ])
273                     [ 0.375    5   ]  [ - 0.625 ]
274 @end group
275 @group
276 (%i4) dlange (inf_norm, b - A . x);
277                        [  1.75   ]
278 (%o4) dlange(inf_norm, [         ]
279                        [ - 0.625 ]
280          [   1    - 2.5 ]         [   1    - 2.5 ]  [  1.75   ]
281        - [              ] . dgesv([              ], [         ]))
282          [ 0.375    5   ]         [ 0.375    5   ]  [ - 0.625 ]
283 @end group
284 @end example
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}.
289 @c ===beg===
290 @c A : matrix ([1, -0.15], [1.82, 2]);
291 @c b : matrix ([3.7, 1, 8], [-2.3, 5, -3.9]);
292 @c x : dgesv (A, b);
293 @c dlange (inf_norm, b - A . x);
294 @c ===end===
295 @example
296 @group
297 (%i1) A : matrix ([1, -0.15], [1.82, 2]);
298                         [  1    - 0.15 ]
299 (%o1)                   [              ]
300                         [ 1.82    2    ]
301 @end group
302 @group
303 (%i2) b : matrix ([3.7, 1, 8], [-2.3, 5, -3.9]);
304                        [  3.7   1    8   ]
305 (%o2)                  [                 ]
306                        [ - 2.3  5  - 3.9 ]
307 @end group
308 @group
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 ]
313 @end group
314 @group
315 (%i4) dlange (inf_norm, b - A . x);
316                        [  3.7   1    8   ]
317 (%o4) dlange(inf_norm, [                 ]
318                        [ - 2.3  5  - 3.9 ]
319    [  1    - 0.15 ]         [  1    - 0.15 ]
320  - [              ] . dgesv([              ], 
321    [ 1.82    2    ]         [ 1.82    2    ]
322 [  3.7   1    8   ]
323 [                 ]))
324 [ - 2.3  5  - 3.9 ]
325 @end group
326 @end example
328 The elements of @var{A} and @var{b} must evaluate to real floating point numbers.
330 @c ===beg===
331 @c A : matrix ([5, -%pi], [1b0, 11/17]);
332 @c b : matrix ([%e], [sin(1)]);
333 @c x : dgesv (A, b);
334 @c dlange (inf_norm, b - A . x);
335 @c ===end===
336 @example
337 @group
338 (%i1) A : matrix ([5, -%pi], [1b0, 11/17]);
339                         [   5    - %pi ]
340                         [              ]
341 (%o1)                   [         11   ]
342                         [ 1.0b0   --   ]
343                         [         17   ]
344 @end group
345 @group
346 (%i2) b : matrix ([%e], [sin(1)]);
347                            [   %e   ]
348 (%o2)                      [        ]
349                            [ sin(1) ]
350 @end group
351 @group
352 (%i3) x : dgesv (A, b);
353                      [   5    - %pi ]
354                      [              ]  [   %e   ]
355 (%o3)          dgesv([         11   ], [        ])
356                      [ 1.0b0   --   ]  [ sin(1) ]
357                      [         17   ]
358 @end group
359 @group
360 (%i4) dlange (inf_norm, b - A . x);
361                        [   %e   ]
362 (%o4) dlange(inf_norm, [        ]
363                        [ sin(1) ]
364           [   5    - %pi ]         [   5    - %pi ]
365           [              ]         [              ]  [   %e   ]
366         - [         11   ] . dgesv([         11   ], [        ]))
367           [ 1.0b0   --   ]         [ 1.0b0   --   ]  [ sin(1) ]
368           [         17   ]         [         17   ]
369 @end group
370 @end example
372 @opencatbox{Categories:}
373 @category{Package lapack}
374 @category{Linear equations}
375 @closecatbox
377 @end deffn
379 @anchor{dgesvd}
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}>>>)
399 such that
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}
404 m4_displaymath(
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>>>)
410 where
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}>>>)
423 that is,
424 m4_mathdot(<<<\mathbf{\Sigma}_{ii} = \sigma_i>>>,
425 <<<@math{@var{Sigma}[i, i]  = @var{sigma}[i]}>>>,
426 <<<{\bf \Sigma}_{ii}  = \sigma_i>>>)
427 The elements
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.
432 The first
433 m4_math(<<<\min(m, n)>>>, <<<min(m, n)>>>)
434 columns of
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}>>>)
444 itself.
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.
456 Example:
458 @c ===beg===
459 @c load ("lapack")$
460 @c fpprintprec : 6;
461 @c M: matrix([1, 2, 3], [3.5, 0.5, 8], [-1, 2, -3], [4, 9, 7]);
462 @c dgesvd (M);
463 @c [sigma, U, VT] : dgesvd (M, true, true);
464 @c m : length (U);
465 @c n : length (VT);
466 @c Sigma:
467 @c   genmatrix(lambda ([i, j], if i=j then sigma[i] else 0),
468 @c             m, n);
469 @c U . Sigma . VT - M;
470 @c transpose (U) . U;
471 @c VT . transpose (VT);
472 @c ===end===
473 @example
474 (%i1) load ("lapack")$
475 @group
476 (%i2) fpprintprec : 6;
477 (%o2)                           6
478 @end group
479 @group
480 (%i3) M: matrix([1, 2, 3], [3.5, 0.5, 8], [-1, 2, -3], [4, 9, 7]);
481                         [  1    2    3  ]
482                         [               ]
483                         [ 3.5  0.5   8  ]
484 (%o3)                   [               ]
485                         [ - 1   2   - 3 ]
486                         [               ]
487                         [  4    9    7  ]
488 @end group
489 @group
490 (%i4) dgesvd (M);
491 (%o4)     [[14.4744, 6.38637, 0.452547], false, false]
492 @end group
493 @group
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 ]
497 [                                                 ]
498 [ - 0.526456   0.672116   - 0.206236   - 0.478091 ]
499 [                                                 ], 
500 [  0.107997   - 0.532278  - 0.0708315  - 0.83666  ]
501 [                                                 ]
502 [ - 0.803287  - 0.514659  - 0.180867    0.239046  ]
503 [ - 0.374486  - 0.538209  - 0.755044 ]
504 [                                    ]
505 [  0.130623   - 0.836799    0.5317   ]]
506 [                                    ]
507 [ - 0.917986   0.100488    0.383672  ]
508 @end group
509 @group
510 (%i6) m : length (U);
511 (%o6)                           4
512 @end group
513 @group
514 (%i7) n : length (VT);
515 (%o7)                           3
516 @end group
517 @group
518 (%i8) Sigma:
519   genmatrix(lambda ([i, j], if i=j then sigma[i] else 0),
520             m, n);
521                  [ 14.4744     0        0     ]
522                  [                            ]
523                  [    0     6.38637     0     ]
524 (%o8)            [                            ]
525                  [    0        0     0.452547 ]
526                  [                            ]
527                  [    0        0        0     ]
528 @end group
529 @group
530 (%i9) U . Sigma . VT - M;
531           [  1.11022e-15        0.0       1.77636e-15 ]
532           [                                           ]
533           [  1.33227e-15    1.66533e-15       0.0     ]
534 (%o9)     [                                           ]
535           [ - 4.44089e-16  - 8.88178e-16  4.44089e-16 ]
536           [                                           ]
537           [  8.88178e-16    1.77636e-15   8.88178e-16 ]
538 @end group
539 @group
540 (%i10) transpose (U) . U;
541        [     1.0      5.55112e-17    2.498e-16     2.77556e-17  ]
542        [                                                        ]
543        [ 5.55112e-17      1.0       5.55112e-17    4.16334e-17  ]
544 (%o10) [                                                        ]
545        [  2.498e-16   5.55112e-17       1.0       - 2.08167e-16 ]
546        [                                                        ]
547        [ 2.77556e-17  4.16334e-17  - 2.08167e-16       1.0      ]
548 @end group
549 @group
550 (%i11) VT . transpose (VT);
551           [      1.0           0.0      - 5.55112e-17 ]
552           [                                           ]
553 (%o11)    [      0.0           1.0       5.55112e-17  ]
554           [                                           ]
555           [ - 5.55112e-17  5.55112e-17       1.0      ]
556 @end group
557 @end example
559 @opencatbox{Categories:}
560 @category{Package lapack}
561 @closecatbox
563 @end deffn
565 @anchor{dlange}
566 @anchor{zlange}
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:
578 @table @code
579 @item max
580 Compute
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.
586 @item one_norm
587 Compute the
588 m4_math(<<<L_1>>>, <<<L[1]>>>)
589 norm of
590 m4_mathcomma(<<<{\bf A}>>>, <<<A>>>)
591 that is, the maximum of the sum of the absolute value of elements in each column.
592 @item inf_norm
593 Compute the
594 m4_math(<<<L_\infty>>>, <<<L[inf]>>>)
595 norm of
596 m4_mathcomma(<<<{\bf A}>>>, <<<A>>>)
597 that is, the maximum of the sum of the absolute value of elements in each row.
598 @item frobenius
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.
602 @end table
604 @opencatbox{Categories:}
605 @category{Package lapack}
606 @closecatbox
608 @end deffn
610 @anchor{dgemm}
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
614 a third matrix.
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}>>>)
625 where
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}>>>)
636 and/or
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:
644 @table @code
645 @item C
646 The matrix
647 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
648 that should be added.  The default is @code{false},
649 which means no matrix is added.
650 @item alpha
651 The product of
652 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
654 m4_math(<<<{\bf B}>>>, <<<@var{B}>>>)
655 is multiplied by this value.  The
656 default is 1.
657 @item beta
658 If a matrix
659 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
660 is given, this value multiplies
661 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
662 before it
663 is added.  The default value is 0, which implies that
664 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
665 is not
666 added, even if
667 m4_math(<<<{\bf C}>>>, <<<@var{C}>>>)
668 is given.  Hence, be sure to specify a non-zero
669 value for
670 m4_mathdot(<<<\beta>>>, <<<@math{beta}>>>)
671 @item transpose_a
672 If @code{true}, the transpose of
673 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
674 is used instead of
675 m4_math(<<<{\bf A}>>>, <<<@var{A}>>>)
676 for the product.  The default is @code{false}.
677 @item transpose_b
678 If @code{true}, the transpose of
679 m4_math(<<<{\bf B}>>>, <<<@var{B}>>>)
680 is used instead of
681 m4_math(<<<{\bf B}>>>, <<<@var{B}>>>)
682 for the product.  The default is @code{false}.
683 @end table
685 @c ===beg===
686 @c load ("lapack")$
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]);
690 @c dgemm(A,B);
691 @c A . B;
692 @c dgemm(A,B,transpose_a=true);
693 @c transpose(A) . B;
694 @c dgemm(A,B,c=C,beta=1);
695 @c A . B + C;
696 @c dgemm(A,B,c=C,beta=1, alpha=-1);
697 @c -A . B + C;
698 @c ===end===
699 @example
700 (%i1) load ("lapack")$
701 @group
702 (%i2) A : matrix([1,2,3],[4,5,6],[7,8,9]);
703                            [ 1  2  3 ]
704                            [         ]
705 (%o2)                      [ 4  5  6 ]
706                            [         ]
707                            [ 7  8  9 ]
708 @end group
709 @group
710 (%i3) B : matrix([-1,-2,-3],[-4,-5,-6],[-7,-8,-9]);
711                         [ - 1  - 2  - 3 ]
712                         [               ]
713 (%o3)                   [ - 4  - 5  - 6 ]
714                         [               ]
715                         [ - 7  - 8  - 9 ]
716 @end group
717 @group
718 (%i4) C : matrix([3,2,1],[6,5,4],[9,8,7]);
719                            [ 3  2  1 ]
720                            [         ]
721 (%o4)                      [ 6  5  4 ]
722                            [         ]
723                            [ 9  8  7 ]
724 @end group
725 @group
726 (%i5) dgemm(A,B);
727                   [ - 30.0   - 36.0   - 42.0  ]
728                   [                           ]
729 (%o5)             [ - 66.0   - 81.0   - 96.0  ]
730                   [                           ]
731                   [ - 102.0  - 126.0  - 150.0 ]
732 @end group
733 @group
734 (%i6) A . B;
735                      [ - 30   - 36   - 42  ]
736                      [                     ]
737 (%o6)                [ - 66   - 81   - 96  ]
738                      [                     ]
739                      [ - 102  - 126  - 150 ]
740 @end group
741 @group
742 (%i7) dgemm(A,B,transpose_a=true);
743                   [ - 66.0  - 78.0   - 90.0  ]
744                   [                          ]
745 (%o7)             [ - 78.0  - 93.0   - 108.0 ]
746                   [                          ]
747                   [ - 90.0  - 108.0  - 126.0 ]
748 @end group
749 @group
750 (%i8) transpose(A) . B;
751                      [ - 66  - 78   - 90  ]
752                      [                    ]
753 (%o8)                [ - 78  - 93   - 108 ]
754                      [                    ]
755                      [ - 90  - 108  - 126 ]
756 @end group
757 @group
758 (%i9) dgemm(A,B,c=C,beta=1);
759                   [ - 27.0  - 34.0   - 41.0  ]
760                   [                          ]
761 (%o9)             [ - 60.0  - 76.0   - 92.0  ]
762                   [                          ]
763                   [ - 93.0  - 118.0  - 143.0 ]
764 @end group
765 @group
766 (%i10) A . B + C;
767                      [ - 27  - 34   - 41  ]
768                      [                    ]
769 (%o10)               [ - 60  - 76   - 92  ]
770                      [                    ]
771                      [ - 93  - 118  - 143 ]
772 @end group
773 @group
774 (%i11) dgemm(A,B,c=C,beta=1, alpha=-1);
775                      [ 33.0   38.0   43.0  ]
776                      [                     ]
777 (%o11)               [ 72.0   86.0   100.0 ]
778                      [                     ]
779                      [ 111.0  134.0  157.0 ]
780 @end group
781 @group
782 (%i12) -A . B + C;
783                         [ 33   38   43  ]
784                         [               ]
785 (%o12)                  [ 72   86   100 ]
786                         [               ]
787                         [ 111  134  157 ]
788 @end group
789 @end example
790 @opencatbox{Categories:}
791 @category{Package lapack}
792 @closecatbox
794 @end deffn
796 @anchor{zgeev}
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}>>>)
802 is complex.
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}
809 @closecatbox
811 @end deffn
813 @anchor{zheev}
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.
833 @c ===beg===
834 @c load("lapack")$
835 @c M: matrix(
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]);
840 @c zheev(M);
841 @c E: zheev(M,true)$
842 @c E[1];
843 @c E[2];
844 @c ===end===
845 @example
846 (%i1) load("lapack")$
847 @group
848 (%i2) M: matrix(
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 ]
854                [                ]         [                  ]
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  ]
858                [                ]         [                  ]
859                [ 9.5 %i - 8.96  ]         [  2.57 - 2.4 %i   ]
860                  [ - 1.72 %i - 1.98 ]         [ - 9.5 %i - 8.96 ]
861                  [                  ]         [                 ]
862                  [  2.25 - 9.51 %i  ]         [  2.4 %i + 2.57  ]
863          Col 3 = [                  ] Col 4 = [                 ]
864                  [      - 4.82      ]         [ 2.04 %i - 3.24  ]
865                  [                  ]         [                 ]
866                  [ - 2.04 %i - 3.24 ]         [      8.44       ]
867 @end group
868 @group
869 (%i3) zheev(M);
870 (%o3) [[- 16.004746472094734, - 6.764970154793324, 
871                    6.6657114535070985, 25.51400517338097], false]
872 @end group
873 (%i4) E: zheev(M,true)$
874 @group
875 (%i5) E[1];
876 (%o5) [- 16.004746472094737, - 6.764970154793325, 
877                            6.665711453507101, 25.514005173380962]
878 @end group
879 @group
880 (%i6) E[2];
881                [  0.26746505331727455 %i + 0.21754535866650165  ]
882                [                                                ]
883                [  0.002696730886619885 %i + 0.6968836773391712  ]
884 (%o6)  Col 1 = [                                                ]
885                [ - 0.6082406376714117 %i - 0.012106142926979313 ]
886                [                                                ]
887                [              0.15930818580950368               ]
888          [  0.26449374706674444 %i + 0.4773693349937472   ]
889          [                                                ]
890          [ - 0.28523890360316206 %i - 0.14143627420116733 ]
891  Col 2 = [                                                ]
892          [  0.2654607680986639 %i + 0.44678181171841735   ]
893          [                                                ]
894          [               0.5750762708542709               ]
895          [  0.28106497673059216 %i - 0.13352639282451817  ]
896          [                                                ]
897          [  0.28663101328695556 %i - 0.4536971347853274   ]
898  Col 3 = [                                                ]
899          [ - 0.29336843237542953 %i - 0.49549724255410565 ]
900          [                                                ]
901          [               0.5325337537576771               ]
902          [ - 0.5737316575503476 %i - 0.39661467994277055 ]
903          [                                               ]
904          [ 0.018265026190214573 %i + 0.35305577043870173 ]
905  Col 4 = [                                               ]
906          [ 0.16737009000854253 %i + 0.01476684746229564  ]
907          [                                               ]
908          [              0.6002632636961784               ]
909 @end group
910 @end example
912 @opencatbox{Categories:}
913 @category{Package lapack}
914 @closecatbox
916 @end deffn