Update docs to match implementation of $build_and_dump_html_index
[maxima.git] / doc / info / es / lapack.texi
blob7ee838a1a75e00986d4fe0134a9bba912c62eab2
1 @c English version 2011-07-15
2 @menu
3 * Introducción a lapack::
4 * Funciones y variables para lapack::
5 @end menu
7 @node Introducción a lapack, Funciones y variables para lapack, lapack, lapack
8 @section Introducción a lapack
10 @code{lapack} es una traducción automática a Common Lisp (con el programa @code{f2c})  
11 de la librería LAPACK escrita en Fortran.
14 @node Funciones y variables para lapack, , Introducción a lapack, lapack
15 @section Funciones y variables para lapack
17 @deffn {Función} dgeev (@var{A})
18 @deffnx {Función} dgeev (@var{A}, @var{right_p}, @var{left_p})
20 Calcula los autovalores y, opcionalmente, también los autovectores de la matriz @var{A}.
21 Todos los elementos de @var{A} deben ser enteros o números decimales en coma flotante.
22 Además, @var{A} debe ser cuadrada (igual número de filas que de columnas) y puede
23 ser o no simétrica.
25 @code{dgeev(@var{A})} calcula sólo los autovalores de @var{A}. 
26 @code{dgeev(@var{A}, @var{right_p}, @var{left_p})} calcula los autovalores de @var{A}
27 y los autovectores por la derecha cuando @math{@var{right_p} = @code{true}}, y
28 los autovectores por la izquierda cuando @math{@var{left_p} = @code{true}}.
30 La función devuelve una lista de tres elementos.
31 El primer elemento es una lista con los autovalores.
32 El segundo elemento es @code{false} o la matriz de autovectores por la derecha.
33 El tercer elemento es @code{false} o la matriz de autovectores por la izquierda.
35 El autovector por la derecha @math{v(j)} (la @math{j}-ésima columna de la matriz de
36 autovectores por la derecha) satisface
38 @math{A . v(j) = lambda(j) . v(j)}
40 donde @math{lambda(j)} es su autovalor asociado.
42 El autovector por la izquierda @math{u(j)} (la @math{j}-ésima columna de la matriz de
43 autovectores por la izquierda) satisface
45 @math{u(j)**H . A = lambda(j) . u(j)**H}
47 donde @math{u(j)**H} denota la transpuesta conjugada de @math{u(j)}.
49 La función de Maxima @code{ctranspose} calcula la transpuesta conjugada.
51 Los autovectores calculados están normalizados para que su norma
52 euclídea valga 1 y su componente mayor tenga su parte
53 imaginaria igual a cero.
55 Ejemplo:
57 @c ===beg===
58 @c load ("lapack")$
59 @c fpprintprec : 6;
60 @c M : matrix ([9.5, 1.75], [3.25, 10.45]);
61 @c dgeev (M);
62 @c [L, v, u] : dgeev (M, true, true);
63 @c D : apply (diag_matrix, L);
64 @c M . v - v . D;
65 @c transpose (u) . M - D . transpose (u);
66 @c ===end===
67 @example
68 (%i1) load ("lapack")$
69 (%i2) fpprintprec : 6;
70 (%o2)                           6
71 (%i3) M : matrix ([9.5, 1.75], [3.25, 10.45]);
72                          [ 9.5   1.75  ]
73 (%o3)                    [             ]
74                          [ 3.25  10.45 ]
75 (%i4) dgeev (M);
76 (%o4)          [[7.54331, 12.4067], false, false]
77 (%i5) [L, v, u] : dgeev (M, true, true);
78                            [ - .666642  - .515792 ]
79 (%o5) [[7.54331, 12.4067], [                      ], 
80                            [  .745378   - .856714 ]
81                                         [ - .856714  - .745378 ]
82                                         [                      ]]
83                                         [  .515792   - .666642 ]
84 (%i6) D : apply (diag_matrix, L);
85                       [ 7.54331     0    ]
86 (%o6)                 [                  ]
87                       [    0     12.4067 ]
88 (%i7) M . v - v . D;
89                 [      0.0       - 8.88178E-16 ]
90 (%o7)           [                              ]
91                 [ - 8.88178E-16       0.0      ]
92 (%i8) transpose (u) . M - D . transpose (u);
93                      [ 0.0  - 4.44089E-16 ]
94 (%o8)                [                    ]
95                      [ 0.0       0.0      ]
96 @end example
97 @end deffn
101 @deffn {Función} dgeqrf (@var{A})
103 Calcula la descomposición QR de la matriz @var{A}.
104 Todos los elementos de @var{A} deben ser enteros o números
105 reales. No es necesario que @var{A} tenga
106 el mismo número de filas que de columnas.
108 La función devuelve una lista con dos elementos; el primero es
109 la matriz cuadrada ortonormal @var{Q}, con el mismo número de
110 filas que @var{A}, y el segundo es la matriz triangular superior @var{R},
111 de iguales dimensiones que @var{A}. El producto @code{@var{Q} . @var{R}}, 
112 siendo "." el operador de la multiplicación matricial, es igual a
113 @var{A}, ignorando errores de redondeo.
115 @c ===beg===
116 @c load ("lapack")$
117 @c fpprintprec : 6;
118 @c M : matrix ([1, -3.2, 8], [-11, 2.7, 5.9]);
119 @c [q, r] : dgeqrf (M);
120 @c q . r - M;
121 @c mat_norm (%);
122 @c ===end===
123 @example
124 (%i1) load ("lapack") $
125 (%i2) fpprintprec : 6 $
126 (%i3) M : matrix ([1, -3.2, 8], [-11, 2.7, 5.9]) $
127 (%i4) [q, r] : dgeqrf (M);
128        [ - .0905357  .995893  ]
129 (%o4) [[                      ], 
130        [  .995893    .0905357 ]
131                                [ - 11.0454   2.97863   5.15148 ]
132                                [                               ]]
133                                [     0      - 2.94241  8.50131 ]
134 (%i5) q . r - M;
135          [ - 7.77156E-16   1.77636E-15   - 8.88178E-16 ]
136 (%o5)    [                                             ]
137          [      0.0       - 1.33227E-15   8.88178E-16  ]
138 (%i6) mat_norm (%, 1);
139 (%o6)                      3.10862E-15
140 @end example
142 @end deffn
147 @deffn {Función} dgesv (@var{A}, @var{b})
149 Calcula la solución @var{x} de la ecuación @math{@var{A} @var{x} = @var{b}},
150 siendo @var{A} una matriz cuadrada y @var{b} otra matriz con el mismo
151 número de filas que @var{A} y un número arbitrario de columnas. Las 
152 dimensiones de la solución @var{x} son las mismas de @var{b}.
154 Los elementos de @var{A} y @var{b} deben ser reducibles a números decimales
155 si se les aplica la función @code{float}, por lo que tales elementos
156 pueden en principio ser de cualquier tipo numérico, constantes numéricas
157 simbólicas o cualesquiera expresiones reducibles a un número decimal.
158 Los elementos de @var{x} son siempre números decimales. Todas las
159 operaciones aritméticas se realizan en coma flotante.
161 @code{dgesv} calcula la solución mediante la descomposición LU de @var{A}.
163 Ejemplos:
165 @code{dgesv} calcula la solución @var{x} de la ecuación 
166 @math{@var{A} @var{x} = @var{b}}.
168 @c ===beg===
169 @c A : matrix ([1, -2.5], [0.375, 5]);
170 @c b : matrix ([1.75], [-0.625]);
171 @c x : dgesv (A, b);
172 @c dlange (inf_norm, b - A . x);
173 @c ===end===
174 @example
175 (%i1) A : matrix ([1, -2.5], [0.375, 5]);
176                                [   1    - 2.5 ]
177 (%o1)                          [              ]
178                                [ 0.375    5   ]
179 (%i2) b : matrix ([1.75], [-0.625]);
180                                   [  1.75   ]
181 (%o2)                             [         ]
182                                   [ - 0.625 ]
183 (%i3) x : dgesv (A, b);
184                             [  1.210526315789474  ]
185 (%o3)                       [                     ]
186                             [ - 0.215789473684211 ]
187 (%i4) dlange (inf_norm, b - A.x);
188 (%o4)                                 0.0
189 @end example
191 @var{b} una matriz con el mismo número de filas que @var{A}
192 y un número arbitrario de columnas. Las dimensiones de 
193 @var{x} son las mismas de @var{b}.
195 @c ===beg===
196 @c A : matrix ([1, -0.15], [1.82, 2]);
197 @c b : matrix ([3.7, 1, 8], [-2.3, 5, -3.9]);
198 @c x : dgesv (A, b);
199 @c dlange (inf_norm, b - A . x);
200 @c ===end===
201 @example
202 (%o0)                                done
203 (%i1) A : matrix ([1, -0.15], [1.82, 2]);
204                                [  1    - 0.15 ]
205 (%o1)                          [              ]
206                                [ 1.82    2    ]
207 (%i2) b : matrix ([3.7, 1, 8], [-2.3, 5, -3.9]);
208                               [  3.7   1    8   ]
209 (%o2)                         [                 ]
210                               [ - 2.3  5  - 3.9 ]
211 (%i3) x : dgesv (A, b);
212         [  3.103827540695117   1.20985481742191    6.781786185657722  ]
213 (%o3)   [                                                             ]
214         [ - 3.974483062032557  1.399032116146062  - 8.121425428948527 ]
215 (%i4) dlange (inf_norm, b - A . x);
216 (%o4)                       1.1102230246251565E-15
217 @end example
219 Los elementos de @var{A} y @var{b} deben ser reducibles a números decimales.
221 @c ===beg===
222 @c A : matrix ([5, -%pi], [1b0, 11/17]);
223 @c b : matrix ([%e], [sin(1)]);
224 @c x : dgesv (A, b);
225 @c dlange (inf_norm, b - A . x);
226 @c ===end===
227 @example
228 (%i1) A : matrix ([5, -%pi], [1b0, 11/17]);
229                                [   5    - %pi ]
230                                [              ]
231 (%o1)                          [         11   ]
232                                [ 1.0b0   --   ]
233                                [         17   ]
234 (%i2) b : matrix ([%e], [sin(1)]);
235                                   [   %e   ]
236 (%o2)                             [        ]
237                                   [ sin(1) ]
238 (%i3) x : dgesv (A, b);
239                              [ 0.690375643155986 ]
240 (%o3)                        [                   ]
241                              [ 0.233510982552952 ]
242 (%i4) dlange (inf_norm, b - A . x);
243 (%o4)                        2.220446049250313E-16
244 @end example
246 @end deffn
251 @deffn {Función} dgesvd (@var{A})
252 @deffnx {Función} dgesvd (@var{A}, @var{left_p}, @var{right_p})
254 Calcula la descomposición singular (SVD, en inglés) de la matriz @var{A},
255 que contiene los valores singulares y, opcionalmente, los vectores singulares por
256 la derecha o por la izquierda. Todos los elementos de @var{A} deben ser enteros o 
257 números decimales en coma flotante. La matriz @var{A} puede ser cuadrada o no 
258 (igual número de filas que de columnas).
260 Sea @math{m} el número de filas y @math{n} el de columnas de @var{A}.
261 La descomposición singular de @var{A} consiste en calcular tres matrices: 
262 @var{U}, @var{Sigma} y @var{V^T}, tales que
264 @c @math{@var{A} = @var{U} . @var{Sigma} . @var{V^T}}
265 @math{@var{A} = @var{U} . @var{Sigma} . @var{V}^T}
267 donde @var{U} es una matriz unitaria @math{m}-por-@math{m},
268 @var{Sigma} es una matriz diagonal @math{m}-por-@math{n} y
269 @var{V^T} es una matriz unitaria @math{n}-por-@math{n}.
271 Sea @math{sigma[i]} un elemento diagonal de @math{Sigma}, esto es,
272 @math{@var{Sigma}[i, i] = @var{sigma}[i]}. Los elementos @math{sigma[i]}
273 se llaman valores singulares de @var{A}, los cuales son reales y no negativos,
274 siendo devueltos por la función @code{dgesvd} en orden descendente. 
276 Las primeras @math{min(m, n)} columnas de @var{U} y @var{V} son los vectores
277 singulares izquierdo y derecho de @var{A}. Nótese que @code{dgesvd}
278 devuelve la transpuesta de @var{V}, no la propia matriz @var{V}.
280 @code{dgesvd(@var{A})} calcula únicamente los valores singulares de @var{A}.
281 @code{dgesvd(@var{A}, @var{left_p}, @var{right_p})} calcula los valores singulares
282 de @var{A} y los vectores sigulares por la izquierda cuando @math{@var{left_p} = @code{true}},
283 y los vectores sigulares por la derecha cuando @math{@var{right_p} = @code{true}}.
285 La función devuelve una lista de tres elementos.
286 El primer elemento es una lista con los valores singulares.
287 El segundo elemento es @code{false} o la matriz de vectores singulares por la izquierda.
288 El tercer elemento es @code{false} o la matriz de vectores singulares por la derecha.
290 Ejemplo:
292 @c ===beg===
293 @c load ("lapack")$
294 @c fpprintprec : 6;
295 @c M: matrix([1, 2, 3], [3.5, 0.5, 8], [-1, 2, -3], [4, 9, 7]);
296 @c dgesvd (M);
297 @c [sigma, U, VT] : dgesvd (M, true, true);
298 @c m : length (U);
299 @c n : length (VT);
300 @c Sigma:
301 @c   genmatrix(lambda ([i, j], if i=j then sigma[i] else 0),
302 @c             m, n);
303 @c U . Sigma . VT - M;
304 @c transpose (U) . U;
305 @c VT . transpose (VT);
306 @c ===end===
307 @example
308 (%i1) load ("lapack")$
309 (%i2) fpprintprec : 6;
310 (%o2)                           6
311 (%i3) M: matrix([1, 2, 3], [3.5, 0.5, 8], [-1, 2, -3], [4, 9, 7]);
312                         [  1    2    3  ]
313                         [               ]
314                         [ 3.5  0.5   8  ]
315 (%o3)                   [               ]
316                         [ - 1   2   - 3 ]
317                         [               ]
318                         [  4    9    7  ]
319 (%i4) dgesvd (M);
320 (%o4)      [[14.4744, 6.38637, .452547], false, false]
321 (%i5) [sigma, U, VT] : dgesvd (M, true, true);
322 (%o5) [[14.4744, 6.38637, .452547], 
323 [ - .256731  .00816168   .959029    - .119523 ]
324 [                                             ]
325 [ - .526456   .672116   - .206236   - .478091 ]
326 [                                             ], 
327 [  .107997   - .532278  - .0708315  - 0.83666 ]
328 [                                             ]
329 [ - .803287  - .514659  - .180867    .239046  ]
330 [ - .374486  - .538209  - .755044 ]
331 [                                 ]
332 [  .130623   - .836799   0.5317   ]]
333 [                                 ]
334 [ - .917986   .100488    .383672  ]
335 (%i6) m : length (U);
336 (%o6)                           4
337 (%i7) n : length (VT);
338 (%o7)                           3
339 (%i8) Sigma:
340         genmatrix(lambda ([i, j], if i=j then sigma[i] else 0),
341                   m, n);
342                   [ 14.4744     0        0    ]
343                   [                           ]
344                   [    0     6.38637     0    ]
345 (%o8)             [                           ]
346                   [    0        0     .452547 ]
347                   [                           ]
348                   [    0        0        0    ]
349 (%i9) U . Sigma . VT - M;
350           [  1.11022E-15        0.0       1.77636E-15 ]
351           [                                           ]
352           [  1.33227E-15    1.66533E-15       0.0     ]
353 (%o9)     [                                           ]
354           [ - 4.44089E-16  - 8.88178E-16  4.44089E-16 ]
355           [                                           ]
356           [  8.88178E-16    1.77636E-15   8.88178E-16 ]
357 (%i10) transpose (U) . U;
358        [     1.0      5.55112E-17    2.498E-16     2.77556E-17  ]
359        [                                                        ]
360        [ 5.55112E-17      1.0       5.55112E-17    4.16334E-17  ]
361 (%o10) [                                                        ]
362        [  2.498E-16   5.55112E-17       1.0       - 2.08167E-16 ]
363        [                                                        ]
364        [ 2.77556E-17  4.16334E-17  - 2.08167E-16       1.0      ]
365 (%i11) VT . transpose (VT);
366           [      1.0           0.0      - 5.55112E-17 ]
367           [                                           ]
368 (%o11)    [      0.0           1.0       5.55112E-17  ]
369           [                                           ]
370           [ - 5.55112E-17  5.55112E-17       1.0      ]
371 @end example
374 @end deffn
376 @deffn {Función} dlange (@var{norm}, @var{A})
377 @deffnx {Función} zlange (@var{norm}, @var{A})
379 Calcula una norma o seudonorma de la matriz @var{A}.
381 @table @code
382 @item max
383 Calcula @math{max(abs(A(i, j)))}, siendo @math{i} y @math{j} números de filas
384 y columnas, respectivamente, de @var{A}.
385 Nótese que esta función no es una norma matricial.
386 @item one_norm
387 Calcula la norma @math{L[1]} de @var{A},
388 esto es, el máximo de la suma de los valores absolutos de los elementos de cada columna.
389 @item inf_norm
390 Calcula la norma @math{L[inf]} de @var{A},
391 esto es, el máximo de la suma de los valores absolutos de los elementos de cada fila.
392 @item frobenius
393 Calcula la norma de Frobenius de @var{A},
394 esto es, la raíz cuadrada de la suma de los cuadrados de los elementos de
395 la matriz.
396 @end table
397 @end deffn
400 @deffn {Función} dgemm (@var{A}, @var{B})
401 @deffnx {Función} dgemm (@var{A}, @var{B}, @var{options})
402 Calcula el producto de dos matrices y, opcionalmente, suma este producto
403 con una tercera matriz.
405 En su forma más simple, @code{dgemm(@var{A}, @var{B})} calcula el
406 producto de las matrices reales @var{A} y @var{B}.
408 En la segunda forma, @code{dgemm} calcula @math{@var{alpha} *
409 @var{A} * @var{B} + @var{beta} * @var{C}}, donde @var{A}, @var{B} y
410 @var{C} son matrices reales de dimensiones apropiadas, siendo @var{alpha} 
411 y @var{beta} números reales. De forma opcional, tanto @var{A} como @var{B}
412 pueden transponerse antes de calcular su producto. Los parámetros adicionales
413 se pueden especificar en cualquier orden, siendo su sintaxis de la forma
414 @code{clave=valor}. Las claves reconocidas son:
416 @table @code
417 @item C
418 La matriz @var{C} que debe ser sumada. El valor por defecto es @code{false},
419 lo que significa que no se sumará ninguna matriz.
420 @item alpha
421 El producto de @var{A} por @var{B} se multiplicará por este vaalor. El valor
422 por defecto es 1.
423 @item beta
424 Si se da la matriz @var{C}, se multiplicará por este valor antes de ser sumada.
425 El valor por defecto es 0, lo que significa que @var{C} no se suma, incluso estando
426 presente. Por lo tanto, téngase cuidado en especificar un valor no nulo para @var{beta}.
427 @item transpose_a
428 Si toma el valor @code{true}, se utilizará la transpuesta de @var{A}, en lugar de la 
429 propia matriz @var{A}, en el producto. El valor por defecto es @code{false}.
430 @item transpose_b
431 Si toma el valor @code{true}, se utilizará la transpuesta de @var{B}, en lugar de la 
432 propia matriz @var{B}, en el producto. El valor por defecto es @code{false}.
433 @end table
435 @c ===beg===
436 @c load ("lapack")$
437 @c A : matrix([1,2,3],[4,5,6],[7,8,9]);
438 @c B : matrix([-1,-2,-3],[-4,-5,-6],[-7,-8,-9]);
439 @c C : matrix([3,2,1],[6,5,4],[9,8,7]);
440 @c dgemm(A,B);
441 @c A . B;
442 @c dgemm(A,B,transpose_a=true);
443 @c transpose(A) . B;
444 @c dgemm(A,B,c=C,beta=1);
445 @c A . B + C;
446 @c dgemm(A,B,c=C,beta=1, alpha=-1);
447 @c -A . B + C
448 @example
449 (%i1) load ("lapack")$
450 (%i2) A : matrix([1,2,3],[4,5,6],[7,8,9]);
451                                   [ 1  2  3 ]
452                                   [         ]
453 (%o2)                             [ 4  5  6 ]
454                                   [         ]
455                                   [ 7  8  9 ]
456 (%i3) B : matrix([-1,-2,-3],[-4,-5,-6],[-7,-8,-9]);
457                                [ - 1  - 2  - 3 ]
458                                [               ]
459 (%o3)                          [ - 4  - 5  - 6 ]
460                                [               ]
461                                [ - 7  - 8  - 9 ]
462 (%i4) C : matrix([3,2,1],[6,5,4],[9,8,7]);
463                                   [ 3  2  1 ]
464                                   [         ]
465 (%o4)                             [ 6  5  4 ]
466                                   [         ]
467                                   [ 9  8  7 ]
468 (%i5) dgemm(A,B);
469                          [ - 30.0   - 36.0   - 42.0  ]
470                          [                           ]
471 (%o5)                    [ - 66.0   - 81.0   - 96.0  ]
472                          [                           ]
473                          [ - 102.0  - 126.0  - 150.0 ]
474 (%i6) A . B;
475                             [ - 30   - 36   - 42  ]
476                             [                     ]
477 (%o6)                       [ - 66   - 81   - 96  ]
478                             [                     ]
479                             [ - 102  - 126  - 150 ]
480 (%i7) dgemm(A,B,transpose_a=true);
481                          [ - 66.0  - 78.0   - 90.0  ]
482                          [                          ]
483 (%o7)                    [ - 78.0  - 93.0   - 108.0 ]
484                          [                          ]
485                          [ - 90.0  - 108.0  - 126.0 ]
486 (%i8) transpose(A) . B;
487                            [ - 66  - 78   - 90  ]
488                            [                    ]
489 (%o8)                      [ - 78  - 93   - 108 ]
490                            [                    ]
491                            [ - 90  - 108  - 126 ]
492 (%i9) dgemm(A,B,c=C,beta=1);
493                          [ - 27.0  - 34.0   - 41.0  ]
494                          [                          ]
495 (%o9)                    [ - 60.0  - 76.0   - 92.0  ]
496                          [                          ]
497                          [ - 93.0  - 118.0  - 143.0 ]
498 (%i10) A . B + C;
499                             [ - 27  - 34   - 41  ]
500                             [                    ]
501 (%o10)                      [ - 60  - 76   - 92  ]
502                             [                    ]
503                             [ - 93  - 118  - 143 ]
504 (%i11) dgemm(A,B,c=C,beta=1, alpha=-1);
505                             [ 33.0   38.0   43.0  ]
506                             [                     ]
507 (%o11)                      [ 72.0   86.0   100.0 ]
508                             [                     ]
509                             [ 111.0  134.0  157.0 ]
510 (%i12) -A . B + C;
511                                [ 33   38   43  ]
512                                [               ]
513 (%o12)                         [ 72   86   100 ]
514                                [               ]
515                                [ 111  134  157 ]
517 @end example
518 @end deffn