source: git/Singular/LIB/linalg.lib @ 0dd77c2

spielwiese
Last change on this file since 0dd77c2 was f0fb366, checked in by Frank Seelisch <seelisch@…>, 14 years ago
fix for non-working example by Bruns (inversion of matrix over the complex numbers) git-svn-id: file:///usr/local/Singular/svn/trunk@12258 2c84dea3-7e68-4137-9b89-c4e89433aadc
  • Property mode set to 100644
File size: 45.6 KB
Line 
1//GMG last modified: 04/25/2000
2//////////////////////////////////////////////////////////////////////////////
3version="$Id$";
4category="Linear Algebra";
5info="
6LIBRARY:  linalg.lib  Algorithmic Linear Algebra
7AUTHORS:  Ivor Saynisch (ivs@math.tu-cottbus.de)
8@*        Mathias Schulze (mschulze@mathematik.uni-kl.de)
9
10PROCEDURES:
11 inverse(A);          matrix, the inverse of A
12 inverse_B(A);        list(matrix Inv,poly p),Inv*A=p*En ( using busadj(A) )
13 inverse_L(A);        list(matrix Inv,poly p),Inv*A=p*En ( using lift )
14 sym_gauss(A);        symmetric gaussian algorithm
15 orthogonalize(A);    Gram-Schmidt orthogonalization
16 diag_test(A);        test whether A can be diagnolized
17 busadj(A);           coefficients of Adj(E*t-A) and coefficients of det(E*t-A)
18 charpoly(A,v);       characteristic polynomial of A ( using busadj(A) )
19 adjoint(A);          adjoint of A ( using busadj(A) )
20 det_B(A);            determinant of A ( using busadj(A) )
21 gaussred(A);         gaussian reduction: P*A=U*S, S a row reduced form of A
22 gaussred_pivot(A);   gaussian reduction: P*A=U*S, uses row pivoting
23 gauss_nf(A);         gaussian normal form of A
24 mat_rk(A);           rank of constant matrix A
25 U_D_O(A);            P*A=U*D*O, P,D,U,O=permutaion,diag,lower-,upper-triang
26 pos_def(A,i);        test symmetric matrix for positive definiteness
27 hessenberg(M);       Hessenberg form of M
28 eigenvals(M);        eigenvalues with multiplicities of M
29 minipoly(M);         minimal polynomial of M
30 spnf(sp);            normal form of spectrum sp
31 spprint(sp);         print spectrum sp
32 jordan(M);           Jordan data of M
33 jordanbasis(M);      Jordan basis and weight filtration of M
34 jordanmatrix(jd);    Jordan matrix with Jordan data jd
35 jordannf(M);         Jordan normal form of M
36";
37
38LIB "matrix.lib";
39LIB "ring.lib";
40LIB "elim.lib";
41LIB "general.lib";
42//////////////////////////////////////////////////////////////////////////////
43// help functions
44//////////////////////////////////////////////////////////////////////////////
45static proc const_mat(matrix A)
46"RETURN:   1 (0) if A is (is not) a constant matrix"
47{
48  int i;
49  int n=ncols(A);
50  def BR=basering;
51  def @R=changeord("dp,c",BR);
52  setring @R;
53  matrix A=fetch(BR,A);
54  for(i=1;i<=n;i=i+1){
55    if(deg(lead(A)[i])>=1){
56      //"input is not a constant matrix";
57      kill @R;
58      setring BR;
59      return(0);
60    }
61  }
62  kill @R;
63  setring BR;
64  return(1);
65}
66//////////////////////////////////////////////////////////////////////////////
67static proc red(matrix A,int i,int j)
68"USAGE:    red(A,i,j);  A = constant matrix
69          reduces column j with respect to A[i,i] and column i
70          reduces row j with respect to A[i,i] and row i
71RETURN:   matrix
72"
73{
74  module m=module(A);
75
76  if(A[i,i]==0){
77    m[i]=m[i]+m[j];
78    m=module(transpose(matrix(m)));
79    m[i]=m[i]+m[j];
80    m=module(transpose(matrix(m)));
81  }
82
83  A=matrix(m);
84  m[j]=m[j]-(A[i,j]/A[i,i])*m[i];
85  m=module(transpose(matrix(m)));
86  m[j]=m[j]-(A[i,j]/A[i,i])*m[i];
87  m=module(transpose(matrix(m)));
88
89  return(matrix(m));
90}
91
92//////////////////////////////////////////////////////////////////////////////
93proc inner_product(vector v1,vector v2)
94"RETURN:   inner product <v1,v2> "
95{
96  int k;
97  if (nrows(v2)>nrows(v1)) { k=nrows(v2); } else { k=nrows(v1); }
98  return ((transpose(matrix(v1,k,1))*matrix(v2,k,1))[1,1]);
99}
100
101/////////////////////////////////////////////////////////////////////////////
102// user functions
103/////////////////////////////////////////////////////////////////////////////
104
105proc inverse(matrix A, list #)
106"USAGE:    inverse(A [,opt]);  A a square matrix, opt integer
107RETURN:
108@format
109          a matrix:
110          - the inverse matrix of A, if A is invertible;
111          - the 1x1 0-matrix if A is not invertible (in the polynomial ring!).
112          There are the following options:
113          - opt=0 or not given: heuristically best option from below
114          - opt=1 : apply std to (transpose(E,A)), ordering (C,dp).
115          - opt=2 : apply interred (transpose(E,A)), ordering (C,dp).
116          - opt=3 : apply lift(A,E), ordering (C,dp).
117@end format
118NOTE:     parameters and minpoly are allowed; opt=2 is only correct for
119          matrices with entries in a field
120SEE ALSO: inverse_B, inverse_L
121EXAMPLE:  example inverse; shows an example
122"
123{
124//--------------------------- initialization and check ------------------------
125   int ii,u,notInvertible,opt;
126   matrix invA;
127   int db = printlevel-voice+3;      //used for comments
128   def R=basering;
129   string mp = string(minpoly);
130   int n = nrows(A);
131   module M = A;
132   intvec v = option(get);           //get options to reset it later
133   if ( ncols(A)!=n )
134   {
135     ERROR("// ** no square matrix");
136   }
137//----------------------- choose heurisitically best option ------------------
138// This may change later, depending on improvements of the implemantation
139// at the monent we use if opt=0 or opt not given:
140// opt = 1 (std) for everything
141// opt = 2 (interred) for nothing, NOTE: interred is ok for constant matricea
142// opt = 3 (lift) for nothing
143// NOTE: interred is ok for constant matrices only (Beispiele am Ende der lib)
144   if(size(#) != 0) {opt = #[1];}
145   if(opt == 0)
146   {
147      if(npars(R) == 0)                       //no parameters
148      {
149         if( size(ideal(A-jet(A,0))) == 0 )   //constant matrix
150         {opt = 1;}
151         else
152         {opt = 1;}
153      }
154      else {opt = 1;}
155   }
156//------------------------- change ring if necessary -------------------------
157   if( ordstr(R) != "C,dp(nvars(R))" )
158   {
159     u=1;
160     def @R=changeord("C,dp",R);
161     setring @R;
162     module M = fetch(R,M);
163     execute("minpoly="+mp+";");
164   }
165//----------------------------- opt=3: use lift ------------------------------
166   if( opt==3 )
167   {
168      module D2;
169      D2 = lift(M,freemodule(n));
170      if (size(ideal(D2))==0)
171      {                                               //catch error in lift
172         dbprint(db,"// ** matrix is not invertible");
173         setring R;
174         if (u==1) { kill @R;}
175         return(invA);
176      }
177   }
178//-------------- opt = 1 resp. opt = 2: use std resp. interred --------------
179   if( opt==1 or opt==2 )
180   {
181      option(redSB);
182      module B = freemodule(n),M;
183      if(opt == 2)
184      {
185         module D = interred(transpose(B));
186         D = transpose(simplify(D,1));
187      }
188      if(opt == 1)
189      {
190         module D = std(transpose(B));
191         D = transpose(simplify(D,1));
192      }
193      module D2 = D[1..n];
194      module D1 = D[n+1..2*n];
195//----------------------- check if matrix is invertible ----------------------
196      for (ii=1; ii<=n; ii++)
197      {
198         if ( D1[ii] != gen(ii) )
199         {
200            notInvertible = 1;
201            break;
202         }
203      }
204   }
205   option(set,v);
206//------------------ return to basering and return result ---------------------
207   if ( u==1 )
208   {
209      setring R;
210      module D2 = fetch(@R,D2);
211      if( opt==1 or opt==2 )
212      { module D1 = fetch(@R,D1);}
213      kill @R;
214   }
215   if( notInvertible == 1 )
216   {
217     // The matrix A seems to be non-invertible.
218     // Note that there are examples, where this is not true but only due to
219     // inexact computations in the field of reals or complex numbers:
220     // ring r = complex, x, dp;
221     // The following matrix has non-zero determinante but seems non-invertible:
222     // matrix A[3][3] = 1,i,i,0,1,2,1,0,1+i;
223     // For this example, inverse_B yields the correct answer.
224     // So, let's use this as a workaround whenever we have this situation:
225     list myList = inverse_B(A);
226     matrix Try = inverse_B(A)[1];
227     if (myList[2] == poly(1)) { return (Try); }
228     else
229     {
230       dbprint(db,"// ** matrix is not invertible");
231       return(invA);
232     }
233   }
234   else { return(matrix(D2)); }     //matrix invertible with inverse D2
235
236}
237example
238{ "EXAMPLE:"; echo = 2;
239  ring r=0,(x,y,z),lp;
240  matrix A[3][3]=
241   1,4,3,
242   1,5,7,
243   0,4,17;
244  print(inverse(A));"";
245  matrix B[3][3]=
246   y+1,  x+y,    y,
247   z,    z+1,    z,
248   y+z+2,x+y+z+2,y+z+1;
249  print(inverse(B));
250  print(B*inverse(B));
251}
252
253//////////////////////////////////////////////////////////////////////////////
254proc sym_gauss(matrix A)
255"USAGE:    sym_gauss(A);  A = symmetric matrix
256RETURN:   matrix, diagonalisation of A with symmetric gauss algorithm
257EXAMPLE:  example sym_gauss; shows an example"
258{
259  int i,j;
260  int n=nrows(A);
261
262  if (ncols(A)!=n){
263    "// ** input is not a square matrix";;
264    return(A);
265  }
266
267  if(!const_mat(A)){
268    "// ** input is not a constant matrix";
269    return(A);
270  }
271
272  if(deg(std(A-transpose(A))[1])!=-1){
273    "// ** input is not a symmetric matrix";
274    return(A);
275  }
276
277  for(i=1; i<n; i++){
278    for(j=i+1; j<=n; j++){
279      if(A[i,j]!=0){ A=red(A,i,j); }
280    }
281  }
282
283  return(A);
284}
285example
286{"EXAMPLE:"; echo = 2;
287  ring r=0,(x),lp;
288  matrix A[2][2]=1,4,4,15;
289  print(A);
290  print(sym_gauss(A));
291}
292
293//////////////////////////////////////////////////////////////////////////////
294proc orthogonalize(matrix A)
295"USAGE:    orthogonalize(A); A = matrix of constants
296RETURN:    matrix, orthogonal basis of the colum space of A
297EXAMPLE:   example orthogonalize; shows an example "
298{
299  int i,j;
300  int n=ncols(A);
301  poly k;
302
303  if(!const_mat(A)){
304    "// ** input is not a constant matrix";
305    matrix B;
306    return(B);
307  }
308
309  module B=module(interred(A));
310
311  for(i=1;i<=n;i=i+1) {
312    for(j=1;j<i;j=j+1) {
313      k=inner_product(B[j],B[j]);
314      if (k==0) { "Error: vector of length zero"; return(matrix(B)); }
315      B[i]=B[i]-(inner_product(B[i],B[j])/k)*B[j];
316    }
317  }
318
319  return(matrix(B));
320}
321example
322{ "EXAMPLE:"; echo = 2;
323  ring r=0,(x),lp;
324  matrix A[4][4]=5,6,12,4,7,3,2,6,12,1,1,2,6,4,2,10;
325  print(A);
326  print(orthogonalize(A));
327}
328
329////////////////////////////////////////////////////////////////////////////
330proc diag_test(matrix A)
331"USAGE:          diag_test(A); A = const square matrix
332RETURN:   int,  1 if A is diagonalizable,@*
333                0 if not@*
334               -1 if no statement is possible, since A does not split.
335NOTE:     The test works only for split matrices, i.e if eigenvalues of A
336          are in the ground field.
337          Does not work with parameters (uses factorize,gcd).
338EXAMPLE:  example diag_test; shows an example"
339{
340  int i,j;
341  int n     = nrows(A);
342  string mp = string(minpoly);
343  string cs = charstr(basering);
344  int np=0;
345
346  if(ncols(A) != n) {
347    "// input is not a square matrix";
348    return(-1);
349  }
350
351   if(!const_mat(A)){
352    "// input is not a constant matrix";
353    return(-1);
354  }
355
356  //Parameterring wegen factorize nicht erlaubt
357  for(i=1;i<size(cs);i=i+1){
358    if(cs[i]==","){np=np+1;} //Anzahl der Parameter
359  }
360  if(np>0){
361        "// rings with parameters not allowed";
362        return(-1);
363  }
364
365  //speichern des aktuellen Rings
366  def BR=basering;
367  //setze R[t]
368  execute("ring rt=("+charstr(basering)+"),(@t,"+varstr(basering)+"),lp;");
369  execute("minpoly="+mp+";");
370  matrix A=imap(BR,A);
371
372  intvec z;
373  intvec s;
374  poly X;         //characteristisches Polynom
375  poly dXdt;      //Ableitung von X nach t
376  ideal g;              //ggT(X,dXdt)
377  poly b;              //Komponente der Busadjunkten-Matrix
378  matrix E[n][n]; //Einheitsmatrix
379
380  E=E+1;
381  A=E*@t-A;
382  X=det(A);
383
384  matrix Xfactors=matrix(factorize(X,1));        //zerfaellt die Matrtix ?
385  int nf=ncols(Xfactors);
386
387  for(i=1;i<=nf;i++){
388    if(lead(Xfactors[1,i])>=@t^2){
389      //" matrix does not split";
390      setring BR;
391      return(-1);
392    }
393  }
394
395  dXdt=diff(X,@t);
396  g=std(ideal(gcd(X,dXdt)));
397
398  //Busadjunkte
399  z=2..n;
400  for(i=1;i<=n;i++){
401    s=2..n;
402    for(j=1;j<=n;j++){
403      b=det(submat(A,z,s));
404
405      if(0!=reduce(b,g)){
406        //" matrix not diagonalizable";
407        setring BR;
408        return(0);
409      }
410
411      s[j]=j;
412    }
413    z[i]=i;
414  }
415
416  //"Die Matrix ist diagonalisierbar";
417  setring BR;
418  return(1);
419}
420example
421{ "EXAMPLE:"; echo = 2;
422  ring r=0,(x),dp;
423  matrix A[4][4]=6,0,0,0,0,0,6,0,0,6,0,0,0,0,0,6;
424  print(A);
425  diag_test(A);
426}
427
428//////////////////////////////////////////////////////////////////////////////
429proc busadj(matrix A)
430"USAGE:   busadj(A);  A = square matrix (of size nxn)
431RETURN:  list L:
432@format
433         L[1] contains the (n+1) coefficients of the characteristic
434              polynomial X of A, i.e.
435              X = L[1][1]+..+L[1][k]*t^(k-1)+..+(L[1][n+1])*t^n
436         L[2] contains the n (nxn)-matrices Hk which are the coefficients of
437              the busadjoint bA = adjoint(E*t-A) of A, i.e.
438              bA = (Hn-1)*t^(n-1)+...+Hk*t^k+...+H0,  ( Hk=L[2][k+1] )
439@end format
440EXAMPLE: example busadj; shows an example"
441{
442  int k;
443  int n    = nrows(A);
444  matrix E = unitmat(n);
445  matrix H[n][n];
446  matrix B[n][n];
447  list bA, X, L;
448  poly a;
449
450  if(ncols(A) != n) {
451    "input is not a square matrix";
452    return(L);
453  }
454
455  bA   = E;
456  X[1] = 1;
457  for(k=1; k<n; k++){
458    B  = A*bA[1];              //bA[1] is the last H
459    a  = -trace(B)/k;
460    H  = B+a*E;
461    bA = insert(bA,H);
462    X  = insert(X,a);
463  }
464  B = A*bA[1];
465  a = -trace(B)/n;
466  X = insert(X,a);
467
468  L = insert(L,bA);
469  L = insert(L,X);
470  return(L);
471}
472example
473{ "EXAMPLE"; echo = 2;
474  ring r = 0,(t,x),lp;
475  matrix A[2][2] = 1,x2,x,x2+3x;
476  print(A);
477  list L = busadj(A);
478  poly X = L[1][1]+L[1][2]*t+L[1][3]*t2; X;
479  matrix bA[2][2] = L[2][1]+L[2][2]*t;
480  print(bA);               //the busadjoint of A;
481  print(bA*(t*unitmat(2)-A));
482}
483
484//////////////////////////////////////////////////////////////////////////////
485proc charpoly(matrix A, list #)
486"USAGE:   charpoly(A[,v]); A square matrix, v string, name of a variable
487RETURN:  poly, the characteristic polynomial det(E*v-A)
488         (default: v=name of last variable)
489NOTE:    A must be independent of the variable v. The computation uses det.
490         If printlevel>0, det(E*v-A) is diplayed recursively.
491EXAMPLE: example charpoly; shows an example"
492{
493  int n = nrows(A);
494  int z = nvars(basering);
495  int i,j;
496  string v;
497  poly X;
498  if(ncols(A) != n)
499  {
500    "// input is not a square matrix";
501    return(X);
502  }
503  //---------------------- test for correct variable -------------------------
504  if( size(#)==0 ){
505    #[1] = varstr(z);
506  }
507  if( typeof(#[1]) == "string") { v = #[1]; }
508  else
509  {
510    "// 2nd argument must be a name of a variable not contained in the matrix";
511    return(X);
512  }
513  j=-1;
514  for(i=1; i<=z; i++)
515  {
516    if(varstr(i)==v){j=i;}
517  }
518  if(j==-1)
519  {
520    "// "+v+" is not a variable in the basering";
521    return(X);
522  }
523  if ( size(select1(module(A),j)) != 0 )
524  {
525    "// matrix must not contain the variable "+v;
526    "// change to a ring with an extra variable, if necessary."
527    return(X);
528  }
529  //-------------------------- compute charpoly ------------------------------
530  X = det(var(j)*unitmat(n)-A);
531  if( printlevel-voice+2 >0) { showrecursive(X,var(j));}
532  return(X);
533}
534example
535{ "EXAMPLE"; echo=2;
536  ring r=0,(x,t),dp;
537  matrix A[3][3]=1,x2,x,x2,6,4,x,4,1;
538  print(A);
539  charpoly(A,"t");
540}
541
542//////////////////////////////////////////////////////////////////////////////
543proc charpoly_B(matrix A, list #)
544"USAGE:   charpoly_B(A[,v]); A square matrix, v string, name of a variable
545RETURN:  poly, the characteristic polynomial det(E*v-A)
546         (default: v=name of last variable)
547NOTE:    A must be constant in the variable v. The computation uses busadj(A).
548EXAMPLE: example charpoly_B; shows an example"
549{
550  int i,j;
551  string s,v;
552  list L;
553  int n     = nrows(A);
554  poly X    = 0;
555  def BR    = basering;
556  string mp = string(minpoly);
557
558  if(ncols(A) != n){
559    "// input is not a square matrix";
560    return(X);
561  }
562
563  //test for correct variable
564  if( size(#)==0 ){
565    #[1] = varstr(nvars(BR));
566  }
567  if( typeof(#[1]) == "string"){
568     v = #[1];
569  }
570  else{
571    "// 2nd argument must be a name of a variable not contained in the matrix";
572    return(X);
573  }
574
575  j=-1;
576  for(i=1; i<=nvars(BR); i++){
577    if(varstr(i)==v){j=i;}
578  }
579  if(j==-1){
580    "// "+v+" is not a variable in the basering";
581    return(X);
582  }
583
584  //var cannot be in A
585  s="Wp(";
586  for( i=1; i<=nvars(BR); i++ ){
587    if(i!=j){ s=s+"0";}
588    else{ s=s+"1";}
589    if( i<nvars(BR)) {s=s+",";}
590  }
591  s=s+")";
592
593  def @R=changeord(s);
594  setring @R;
595  execute("minpoly="+mp+";");
596  matrix A = imap(BR,A);
597  for(i=1; i<=n; i++){
598    if(deg(lead(A)[i])>=1){
599      "// matrix must not contain the variable "+v;
600      kill @R;
601      setring BR;
602      return(X);
603    }
604  }
605
606  //get coefficients and build the char. poly
607  kill @R;
608  setring BR;
609  L = busadj(A);
610  for(i=1; i<=n+1; i++){
611    execute("X=X+L[1][i]*"+v+"^"+string(i-1)+";");
612  }
613
614  return(X);
615}
616example
617{ "EXAMPLE"; echo=2;
618  ring r=0,(x,t),dp;
619  matrix A[3][3]=1,x2,x,x2,6,4,x,4,1;
620  print(A);
621  charpoly_B(A,"t");
622}
623
624//////////////////////////////////////////////////////////////////////////////
625proc adjoint(matrix A)
626"USAGE:    adjoint(A);  A = square matrix
627RETURN:   adjoint matrix of A, i.e. Adj*A=det(A)*E
628NOTE:     computation uses busadj(A)
629EXAMPLE:  example adjoint; shows an example"
630{
631  int n=nrows(A);
632  matrix Adj[n][n];
633  list L;
634
635  if(ncols(A) != n) {
636    "// input is not a square matrix";
637    return(Adj);
638  }
639
640  L  = busadj(A);
641  Adj= (-1)^(n-1)*L[2][1];
642  return(Adj);
643
644}
645example
646{ "EXAMPLE"; echo=2;
647  ring r=0,(t,x),lp;
648  matrix A[2][2]=1,x2,x,x2+3x;
649  print(A);
650  matrix Adj[2][2]=adjoint(A);
651  print(Adj);                    //Adj*A=det(A)*E
652  print(Adj*A);
653}
654
655//////////////////////////////////////////////////////////////////////////////
656proc inverse_B(matrix A)
657"USAGE:    inverse_B(A);  A = square matrix
658RETURN:   list Inv with
659          - Inv[1] = matrix I and
660          - Inv[2] = poly p
661          such that I*A = unitmat(n)*p;
662NOTE:     p=1 if 1/det(A) is computable and  p=det(A) if not;
663          the computation uses busadj.
664SEE ALSO: inverse, inverse_L
665EXAMPLE:  example inverse_B; shows an example"
666{
667  int i;
668  int n=nrows(A);
669  matrix I[n][n];
670  poly factor;
671  list L;
672  list Inv;
673
674  if(ncols(A) != n) {
675    "input is not a square matrix";
676    return(I);
677  }
678
679  L=busadj(A);
680  I=module(-L[2][1]);        //+-Adj(A)
681
682  if(reduce(1,std(L[1][1]))==0){
683    I=I*lift(L[1][1],1)[1][1];
684    factor=1;
685  }
686  else{ factor=L[1][1];}     //=+-det(A) or 1
687  Inv=insert(Inv,factor);
688  Inv=insert(Inv,matrix(I));
689
690  return(Inv);
691}
692example
693{ "EXAMPLE"; echo=2;
694  ring r=0,(x,y),lp;
695  matrix A[3][3]=x,y,1,1,x2,y,x,6,0;
696  print(A);
697  list Inv=inverse_B(A);
698  print(Inv[1]);
699  print(Inv[2]);
700  print(Inv[1]*A);
701}
702
703//////////////////////////////////////////////////////////////////////////////
704proc det_B(matrix A)
705"USAGE:     det_B(A);  A any matrix
706RETURN:    returns the determinant of A
707NOTE:      the computation uses the busadj algorithm
708EXAMPLE:   example det_B; shows an example"
709{
710  int n=nrows(A);
711  list L;
712
713  if(ncols(A) != n){ return(0);}
714
715  L=busadj(A);
716  return((-1)^n*L[1][1]);
717}
718example
719{ "EXAMPLE"; echo=2;
720  ring r=0,(x),dp;
721  matrix A[10][10]=random(2,10,10)+unitmat(10)*x;
722  print(A);
723  det_B(A);
724}
725
726//////////////////////////////////////////////////////////////////////////////
727proc inverse_L(matrix A)
728"USAGE:     inverse_L(A);  A = square matrix
729RETURN:    list Inv representing a left inverse of A, i.e
730           - Inv[1] = matrix I and
731           - Inv[2] = poly p
732           such that I*A = unitmat(n)*p;
733NOTE:      p=1 if 1/det(A) is computable and p=det(A) if not;
734           the computation computes first det(A) and then uses lift
735SEE ALSO:  inverse, inverse_B
736EXAMPLE:   example inverse_L; shows an example"
737{
738  int n=nrows(A);
739  matrix I;
740  matrix E[n][n]=unitmat(n);
741  poly factor;
742  poly d=1;
743  list Inv;
744
745  if (ncols(A)!=n){
746    "// input is not a square matrix";
747    return(I);
748  }
749
750  d=det(A);
751  if(d==0){
752    "// matrix is not invertible";
753    return(Inv);
754  }
755
756  // test if 1/det(A) exists
757  if(reduce(1,std(d))!=0){ E=E*d;}
758
759  I=lift(A,E);
760  if(I==unitmat(n)-unitmat(n)){ //catch error in lift
761    "// matrix is not invertible";
762    return(Inv);
763  }
764
765  factor=d;      //=det(A) or 1
766  Inv=insert(Inv,factor);
767  Inv=insert(Inv,I);
768
769  return(Inv);
770}
771example
772{ "EXAMPLE"; echo=2;
773  ring r=0,(x,y),lp;
774  matrix A[3][3]=x,y,1,1,x2,y,x,6,0;
775  print(A);
776  list Inv=inverse_L(A);
777  print(Inv[1]);
778  print(Inv[2]);
779  print(Inv[1]*A);
780}
781
782//////////////////////////////////////////////////////////////////////////////
783proc gaussred(matrix A)
784"USAGE:   gaussred(A);   A any constant matrix
785RETURN:  list Z:  Z[1]=P , Z[2]=U , Z[3]=S , Z[4]=rank(A)
786         gives a row reduced matrix S, a permutation matrix P and a
787         normalized lower triangular matrix U, with P*A=U*S
788NOTE:    This procedure is designed for teaching purposes mainly.
789         The straight forward implementation in the interpreted library
790         is not very efficient (no standard basis computation).
791EXAMPLE: example gaussred; shows an example"
792{
793  int i,j,l,k,jp,rang;
794  poly c,pivo;
795  list Z;
796  int n = nrows(A);
797  int m = ncols(A);
798  int mr= n; //max. rang
799  matrix P[n][n] = unitmat(n);
800  matrix U[n][n] = P;
801
802  if(!const_mat(A)){
803    "// input is not a constant matrix";
804    return(Z);
805  }
806
807  if(n>m){mr=m;} //max. rang
808
809  for(i=1;i<=mr;i=i+1){
810    if((i+k)>m){break};
811
812    //Test: Diagonalelement=0
813    if(A[i,i+k]==0){
814      jp=i;pivo=0;
815      for(j=i+1;j<=n;j=j+1){
816        c=absValue(A[j,i+k]);
817        if(pivo<c){ pivo=c;jp=j;}
818      }
819      if(jp != i){       //Zeilentausch
820        for(j=1;j<=m;j=j+1){ //Zeilentausch in A (und U) (i-te mit jp-ter)
821          c=A[i,j];
822          A[i,j]=A[jp,j];
823          A[jp,j]=c;
824        }
825        for(j=1;j<=n;j=j+1){ //Zeilentausch in P
826          c=P[i,j];
827          P[i,j]=P[jp,j];
828          P[jp,j]=c;
829        }
830      }
831      if(pivo==0){k++;continue;} //eine von selbst auftauchende Stufe !
832    }                          //i sollte im naechsten Lauf nicht erhoeht sein
833
834    //Eliminationsschritt
835    for(j=i+1;j<=n;j=j+1){
836      c=A[j,i+k]/A[i,i+k];
837      for(l=i+k+1;l<=m;l=l+1){
838        A[j,l]=A[j,l]-A[i,l]*c;
839      }
840      A[j,i+k]=0;  // nur wichtig falls k>0 ist
841      A[j,i]=c;    // bildet U
842    }
843  rang=i;
844  }
845
846  for(i=1;i<=mr;i=i+1){
847    for(j=i+1;j<=n;j=j+1){
848      U[j,i]=A[j,i];
849      A[j,i]=0;
850    }
851  }
852
853  Z=insert(Z,rang);
854  Z=insert(Z,A);
855  Z=insert(Z,U);
856  Z=insert(Z,P);
857
858  return(Z);
859}
860example
861{ "EXAMPLE";echo=2;
862  ring r=0,(x),dp;
863  matrix A[5][4]=1,3,-1,4,2,5,-1,3,1,3,-1,4,0,4,-3,1,-3,1,-5,-2;
864  print(A);
865  list Z=gaussred(A);   //construct P,U,S s.t. P*A=U*S
866  print(Z[1]);          //P
867  print(Z[2]);          //U
868  print(Z[3]);          //S
869  print(Z[4]);          //rank
870  print(Z[1]*A);        //P*A
871  print(Z[2]*Z[3]);     //U*S
872}
873
874//////////////////////////////////////////////////////////////////////////////
875proc gaussred_pivot(matrix A)
876"USAGE:     gaussred_pivot(A);   A any constant matrix
877RETURN:    list Z:  Z[1]=P , Z[2]=U , Z[3]=S , Z[4]=rank(A)
878           gives a row reduced matrix S, a permutation matrix P and a
879           normalized lower triangular matrix U, with P*A=U*S
880NOTE:      with row pivoting
881EXAMPLE:   example gaussred_pivot; shows an example"
882{
883  int i,j,l,k,jp,rang;
884  poly c,pivo;
885  list Z;
886  int n=nrows(A);
887  int m=ncols(A);
888  int mr=n; //max. rang
889  matrix P[n][n]=unitmat(n);
890  matrix U[n][n]=P;
891
892  if(!const_mat(A)){
893    "// input is not a constant matrix";
894    return(Z);
895  }
896
897  if(n>m){mr=m;} //max. rang
898
899  for(i=1;i<=mr;i=i+1){
900    if((i+k)>m){break};
901
902    //Pivotisierung
903    pivo=absValue(A[i,i+k]);jp=i;
904    for(j=i+1;j<=n;j=j+1){
905      c=absValue(A[j,i+k]);
906      if(pivo<c){ pivo=c;jp=j;}
907    }
908    if(jp != i){ //Zeilentausch
909      for(j=1;j<=m;j=j+1){ //Zeilentausch in A (und U) (i-te mit jp-ter)
910        c=A[i,j];
911        A[i,j]=A[jp,j];
912        A[jp,j]=c;
913      }
914      for(j=1;j<=n;j=j+1){ //Zeilentausch in P
915        c=P[i,j];
916        P[i,j]=P[jp,j];
917        P[jp,j]=c;
918      }
919    }
920    if(pivo==0){k++;continue;} //eine von selbst auftauchende Stufe !
921                               //i sollte im naechsten Lauf nicht erhoeht sein
922    //Eliminationsschritt
923    for(j=i+1;j<=n;j=j+1){
924      c=A[j,i+k]/A[i,i+k];
925      for(l=i+k+1;l<=m;l=l+1){
926        A[j,l]=A[j,l]-A[i,l]*c;
927      }
928      A[j,i+k]=0;  // nur wichtig falls k>0 ist
929      A[j,i]=c;    // bildet U
930    }
931  rang=i;
932  }
933
934  for(i=1;i<=mr;i=i+1){
935    for(j=i+1;j<=n;j=j+1){
936      U[j,i]=A[j,i];
937      A[j,i]=0;
938    }
939  }
940
941  Z=insert(Z,rang);
942  Z=insert(Z,A);
943  Z=insert(Z,U);
944  Z=insert(Z,P);
945
946  return(Z);
947}
948example
949{ "EXAMPLE";echo=2;
950  ring r=0,(x),dp;
951  matrix A[5][4] = 1, 3,-1,4,
952                   2, 5,-1,3,
953                   1, 3,-1,4,
954                   0, 4,-3,1,
955                  -3,1,-5,-2;
956  list Z=gaussred_pivot(A);  //construct P,U,S s.t. P*A=U*S
957  print(Z[1]);               //P
958  print(Z[2]);               //U
959  print(Z[3]);               //S
960  print(Z[4]);               //rank
961  print(Z[1]*A);             //P*A
962  print(Z[2]*Z[3]);          //U*S
963}
964
965//////////////////////////////////////////////////////////////////////////////
966proc gauss_nf(matrix A)
967"USAGE:    gauss_nf(A); A any constant matrix
968RETURN:   matrix; gauss normal form of A (uses gaussred)
969EXAMPLE:  example gauss_nf; shows an example"
970{
971  list Z;
972  if(!const_mat(A)){
973    "// input is not a constant matrix";
974    return(A);
975  }
976  Z = gaussred(A);
977  return(Z[3]);
978}
979example
980{ "EXAMPLE";echo=2;
981  ring r = 0,(x),dp;
982  matrix A[4][4] = 1,4,4,7,2,5,5,4,4,1,1,3,0,2,2,7;
983  print(gauss_nf(A));
984}
985
986//////////////////////////////////////////////////////////////////////////////
987proc mat_rk(matrix A)
988"USAGE:    mat_rk(A); A any constant matrix
989RETURN:   int, rank of A
990EXAMPLE:  example mat_rk; shows an example"
991{
992  list Z;
993  if(!const_mat(A)){
994    "// input is not a constant matrix";
995    return(-1);
996  }
997  Z = gaussred(A);
998  return(Z[4]);
999}
1000example
1001{ "EXAMPLE";echo=2;
1002  ring r = 0,(x),dp;
1003  matrix A[4][4] = 1,4,4,7,2,5,5,4,4,1,1,3,0,2,2,7;
1004  mat_rk(A);
1005}
1006
1007//////////////////////////////////////////////////////////////////////////////
1008proc U_D_O(matrix A)
1009"USAGE:     U_D_O(A);   constant invertible matrix A
1010RETURN:    list Z:  Z[1]=P , Z[2]=U , Z[3]=D , Z[4]=O
1011           gives a permutation matrix P,
1012           a normalized lower triangular matrix U ,
1013           a diagonal matrix D, and
1014           a normalized upper triangular matrix O
1015           with P*A=U*D*O
1016NOTE:      Z[1]=-1 means that A is not regular (proc uses gaussred)
1017EXAMPLE:   example U_D_O; shows an example"
1018{
1019  int i,j;
1020  list Z,L;
1021  int n=nrows(A);
1022  matrix O[n][n]=unitmat(n);
1023  matrix D[n][n];
1024
1025  if (ncols(A)!=n){
1026    "// input is not a square matrix";
1027    return(Z);
1028  }
1029  if(!const_mat(A)){
1030    "// input is not a constant matrix";
1031    return(Z);
1032  }
1033
1034  L=gaussred(A);
1035
1036  if(L[4]!=n){
1037    "// input is not an invertible matrix";
1038    Z=insert(Z,-1);  //hint for calling procedures
1039    return(Z);
1040  }
1041
1042  D=L[3];
1043
1044  for(i=1; i<=n; i++){
1045    for(j=i+1; j<=n; j++){
1046      O[i,j] = D[i,j]/D[i,i];
1047      D[i,j] = 0;
1048    }
1049  }
1050
1051  Z=insert(Z,O);
1052  Z=insert(Z,D);
1053  Z=insert(Z,L[2]);
1054  Z=insert(Z,L[1]);
1055  return(Z);
1056}
1057example
1058{ "EXAMPLE";echo=2;
1059  ring r = 0,(x),dp;
1060  matrix A[5][5] = 10, 4,  0, -9,  8,
1061                   -3, 6, -6, -4,  9,
1062                    0, 3, -1, -9, -8,
1063                   -4,-2, -6, -10,10,
1064                   -9, 5, -1, -6,  5;
1065  list Z = U_D_O(A);              //construct P,U,D,O s.t. P*A=U*D*O
1066  print(Z[1]);                    //P
1067  print(Z[2]);                    //U
1068  print(Z[3]);                    //D
1069  print(Z[4]);                    //O
1070  print(Z[1]*A);                  //P*A
1071  print(Z[2]*Z[3]*Z[4]);          //U*D*O
1072}
1073
1074//////////////////////////////////////////////////////////////////////////////
1075proc pos_def(matrix A)
1076"USAGE:     pos_def(A); A = constant, symmetric square matrix
1077RETURN:    int:
1078           1  if A is positive definit ,
1079           0  if not,
1080           -1 if unknown
1081EXAMPLE:   example pos_def; shows an example"
1082{
1083  int j;
1084  list Z;
1085  int n = nrows(A);
1086  matrix H[n][n];
1087
1088  if (ncols(A)!=n){
1089    "// input is not a square matrix";
1090    return(0);
1091  }
1092  if(!const_mat(A)){
1093    "// input is not a constant matrix";
1094    return(-1);
1095  }
1096  if(deg(std(A-transpose(A))[1])!=-1){
1097    "// input is not a hermitian (symmetric) matrix";
1098    return(-1);
1099  }
1100
1101  Z=U_D_O(A);
1102
1103  if(Z[1]==-1){
1104    return(0);
1105  }  //A not regular, therefore not pos. definit
1106
1107  H=Z[1];
1108  //es fand Zeilentausch statt: also nicht positiv definit
1109  if(deg(std(H-unitmat(n))[1])!=-1){
1110    return(0);
1111  }
1112
1113  H=Z[3];
1114
1115  for(j=1;j<=n;j=j+1){
1116    if(H[j,j]<=0){
1117      return(0);
1118    } //eigenvalue<=0, not pos.definit
1119  }
1120
1121  return(1); //positiv definit;
1122}
1123example
1124{ "EXAMPLE"; echo=2;
1125  ring r = 0,(x),dp;
1126  matrix A[5][5] = 20,  4,  0, -9,   8,
1127                    4, 12, -6, -4,   9,
1128                    0, -6, -2, -9,  -8,
1129                   -9, -4, -9, -20, 10,
1130                    8,  9, -8,  10, 10;
1131  pos_def(A);
1132  matrix B[3][3] =  3,  2,  0,
1133                    2, 12,  4,
1134                    0,  4,  2;
1135  pos_def(B);
1136}
1137
1138//////////////////////////////////////////////////////////////////////////////
1139proc linsolve(matrix A, matrix b)
1140"USAGE:     linsolve(A,b); A a constant nxm-matrix, b a constant nx1-matrix
1141RETURN:    a 1xm matrix X, solution of inhomogeneous linear system A*X = b
1142           return the 0-matrix if system is not solvable
1143NOTE:      uses gaussred
1144EXAMPLE:   example linsolve; shows an example"
1145{
1146  int i,j,k,rc,r;
1147  poly c;
1148  list Z;
1149  int n  = nrows(A);
1150  int m  = ncols(A);
1151  int n_b= nrows(b);
1152  matrix Ab[n][m+1];
1153  matrix X[m][1];
1154
1155  if(ncols(b)!=1){
1156    "// right hand side b is not a nx1 matrix";
1157    return(X);
1158  }
1159
1160  if(!const_mat(A)){
1161    "// input hand is not a constant matrix";
1162    return(X);
1163  }
1164
1165  if(n_b>n){
1166    for(i=n; i<=n_b; i++){
1167      if(b[i,1]!=0){
1168        "// right hand side b not in Image(A)";
1169        return X;
1170      }
1171    }
1172  }
1173
1174  if(n_b<n){
1175    matrix copy[n_b][1]=b;
1176    matrix b[n][1]=0;
1177    for(i=1;i<=n_b;i=i+1){
1178      b[i,1]=copy[i,1];
1179    }
1180  }
1181
1182  r=mat_rk(A);
1183
1184  //1. b constant vector
1185  if(const_mat(b)){
1186    //extend A with b
1187    for(i=1; i<=n; i++){
1188      for(j=1; j<=m; j++){
1189        Ab[i,j]=A[i,j];
1190      }
1191      Ab[i,m+1]=b[i,1];
1192    }
1193
1194    //Gauss reduction
1195    Z  = gaussred(Ab);
1196    Ab = Z[3];  //normal form
1197    rc = Z[4];  //rank(Ab)
1198    //print(Ab);
1199
1200    if(r<rc){
1201        "// no solution";
1202        return(X);
1203    }
1204    k=m;
1205    for(i=r;i>=1;i=i-1){
1206
1207      j=1;
1208      while(Ab[i,j]==0){j=j+1;}// suche Ecke
1209
1210      for(;k>j;k=k-1){ X[k]=0;}//springe zur Ecke
1211
1212
1213      c=Ab[i,m+1]; //i-te Komponene von b
1214      for(j=m;j>k;j=j-1){
1215        c=c-X[j,1]*Ab[i,j];
1216      }
1217      if(Ab[i,k]==0){
1218        X[k,1]=1; //willkuerlich
1219      }
1220      else{
1221          X[k,1]=c/Ab[i,k];
1222      }
1223      k=k-1;
1224      if(k==0){break;}
1225    }
1226
1227
1228  }//endif (const b)
1229  else{  //b not constant
1230    "// !not implemented!";
1231
1232  }
1233
1234  return(X);
1235}
1236example
1237{ "EXAMPLE";echo=2;
1238  ring r=0,(x),dp;
1239  matrix A[3][2] = -4,-6,
1240                    2, 3,
1241                   -5, 7;
1242  matrix b[3][1] = 10,
1243                   -5,
1244                    2;
1245  matrix X = linsolve(A,b);
1246  print(X);
1247  print(A*X);
1248}
1249//////////////////////////////////////////////////////////////////////////////
1250
1251///////////////////////////////////////////////////////////////////////////////
1252//    PROCEDURES for Jordan normal form
1253//    AUTHOR:  Mathias Schulze, email: mschulze@mathematik.uni-kl.de
1254///////////////////////////////////////////////////////////////////////////////
1255
1256static proc rowcolswap(matrix M,int i,int j)
1257{
1258  if(i==j)
1259  {
1260    return(M);
1261  }
1262  poly p;
1263  for(int k=1;k<=nrows(M);k++)
1264  {
1265    p=M[i,k];
1266    M[i,k]=M[j,k];
1267    M[j,k]=p;
1268  }
1269  for(k=1;k<=ncols(M);k++)
1270  {
1271    p=M[k,i];
1272    M[k,i]=M[k,j];
1273    M[k,j]=p;
1274  }
1275  return(M);
1276}
1277//////////////////////////////////////////////////////////////////////////////
1278
1279static proc rowelim(matrix M,int i,int j,int k)
1280{
1281  if(jet(M[i,k],0)==0||jet(M[j,k],0)==0)
1282  {
1283    return(M);
1284  }
1285  number n=number(jet(M[i,k],0))/number(jet(M[j,k],0));
1286  for(int l=1;l<=ncols(M);l++)
1287  {
1288    M[i,l]=M[i,l]-n*M[j,l];
1289  }
1290  for(l=1;l<=nrows(M);l++)
1291  {
1292    M[l,j]=M[l,j]+n*M[l,i];
1293  }
1294  return(M);
1295}
1296///////////////////////////////////////////////////////////////////////////////
1297
1298static proc colelim(matrix M,int i,int j,int k)
1299{
1300  if(jet(M[k,i],0)==0||jet(M[k,j],0)==0)
1301  {
1302    return(M);
1303  }
1304  number n=number(jet(M[k,i],0))/number(jet(M[k,j],0));
1305  for(int l=1;l<=nrows(M);l++)
1306  {
1307    M[l,i]=M[l,i]-n*M[l,j];
1308  }
1309  for(l=1;l<=ncols(M);l++)
1310  {
1311    M[j,l]=M[j,l]+n*M[i,l];
1312  }
1313  return(M);
1314}
1315///////////////////////////////////////////////////////////////////////////////
1316
1317proc hessenberg(matrix M)
1318"USAGE:   hessenberg(M); matrix M
1319ASSUME:  M constant square matrix
1320RETURN:  matrix H;  Hessenberg form of M
1321EXAMPLE: example hessenberg; shows examples
1322"
1323{
1324  if(system("with","eigenval"))
1325  {
1326    return(system("hessenberg",M));
1327  }
1328
1329  int n=ncols(M);
1330  int i,j;
1331  for(int k=1;k<n-1;k++)
1332  {
1333    j=k+1;
1334    while(j<n&&jet(M[j,k],0)==0)
1335    {
1336      j++;
1337    }
1338    if(jet(M[j,k],0)!=0)
1339    {
1340      M=rowcolswap(M,j,k+1);
1341      for(i=j+1;i<=n;i++)
1342      {
1343        M=rowelim(M,i,k+1,k);
1344      }
1345    }
1346  }
1347  return(M);
1348}
1349example
1350{ "EXAMPLE:"; echo=2;
1351  ring R=0,x,dp;
1352  matrix M[3][3]=3,2,1,0,2,1,0,0,3;
1353  print(M);
1354  print(hessenberg(M));
1355}
1356///////////////////////////////////////////////////////////////////////////////
1357
1358proc eigenvals(matrix M)
1359"USAGE:   eigenvals(M); matrix M
1360ASSUME:  eigenvalues of M in basefield
1361RETURN:
1362@format
1363list l;
1364  ideal l[1];
1365    number l[1][i];  i-th eigenvalue of M
1366  intvec l[2];
1367    int l[2][i];  multiplicity of i-th eigenvalue of M
1368@end format
1369EXAMPLE: example eigenvals; shows examples
1370"
1371{
1372  if(system("with","eigenval"))
1373  {
1374    return(system("eigenvals",jet(M,0)));
1375  }
1376
1377  M=jet(hessenberg(M),0);
1378  int n=ncols(M);
1379  int k;
1380  ideal e;
1381  intvec m;
1382  number e0;
1383  intvec v;
1384  list l;
1385  int i,j;
1386  j=1;
1387  while(j<=n)
1388  {
1389    v=j;
1390    j++;
1391    if(j<=n)
1392    {
1393      while(j<n&&M[j,j-1]!=0)
1394      {
1395        v=v,j;
1396        j++;
1397      }
1398      if(M[j,j-1]!=0)
1399      {
1400        v=v,j;
1401        j++;
1402      }
1403    }
1404    if(size(v)==1)
1405    {
1406      k++;
1407      e[k]=M[v,v];
1408      m[k]=1;
1409    }
1410    else
1411    {
1412      l=factorize(det(submat(M,v,v)-var(1)));
1413      for(i=size(l[1]);i>=1;i--)
1414      {
1415        e0=number(jet(l[1][i]/var(1),0));
1416        if(e0!=0)
1417        {
1418          k++;
1419          e[k]=(e0*var(1)-l[1][i])/e0;
1420          m[k]=l[2][i];
1421        }
1422      }
1423    }
1424  }
1425  return(spnf(list(e,m)));
1426}
1427example
1428{ "EXAMPLE:"; echo=2;
1429  ring R=0,x,dp;
1430  matrix M[3][3]=3,2,1,0,2,1,0,0,3;
1431  print(M);
1432  eigenvals(M);
1433}
1434///////////////////////////////////////////////////////////////////////////////
1435
1436proc minipoly(matrix M,list #)
1437"USAGE:   minipoly(M); matrix M
1438ASSUME:  eigenvalues of M in basefield
1439RETURN:
1440@format
1441list l;  minimal polynomial of M
1442  ideal l[1];
1443    number l[1][i];  i-th root of minimal polynomial of M
1444  intvec l[2];
1445    int l[2][i];  multiplicity of i-th root of minimal polynomial of M
1446@end format
1447EXAMPLE: example minipoly; shows examples
1448"
1449{
1450  if(nrows(M)==0)
1451  {
1452    ERROR("non empty expected");
1453  }
1454  if(ncols(M)!=nrows(M))
1455  {
1456    ERROR("square matrix expected");
1457  }
1458
1459  M=jet(M,0);
1460
1461  if(size(#)==0)
1462  {
1463    #=eigenvals(M);
1464  }
1465  def e0,m0=#[1..2];
1466
1467  intvec m1;
1468  matrix N0,N1;
1469  for(int i=1;i<=ncols(e0);i++)
1470  {
1471    m1[i]=1;
1472    N0=M-e0[i];
1473    N1=N0;
1474    while(size(syz(N1))<m0[i])
1475    {
1476      m1[i]=m1[i]+1;
1477      N1=N1*N0;
1478    }
1479  }
1480
1481  return(list(e0,m1));
1482}
1483example
1484{ "EXAMPLE:"; echo=2;
1485  ring R=0,x,dp;
1486  matrix M[3][3]=3,2,1,0,2,1,0,0,3;
1487  print(M);
1488  minipoly(M);
1489}
1490///////////////////////////////////////////////////////////////////////////////
1491
1492proc spnf(list #)
1493"USAGE:   spnf(list(a[,m])); ideal a, intvec m
1494ASSUME:  ncols(a)==size(m)
1495RETURN:  list l:
1496            l[1] an ideal, the generators of a; sorted and with multiple entries displayed only once@*
1497            l[2] and intvec, l[2][i] provides the multiplicity of l[1][i]
1498EXAMPLE: example spnf; shows examples
1499"
1500{
1501  list sp=#;
1502  ideal a=sp[1];
1503  int n=ncols(a);
1504  intvec m;
1505  list V;
1506  module v;
1507  int i,j;
1508  for(i=2;i<=size(sp);i++)
1509  {
1510    if(typeof(sp[i])=="intvec")
1511    {
1512      m=sp[i];
1513    }
1514    if(typeof(sp[i])=="module")
1515    {
1516      v=sp[i];
1517      for(j=n;j>=1;j--)
1518      {
1519        V[j]=module(v[j]);
1520      }
1521    }
1522    if(typeof(sp[i])=="list")
1523    {
1524      V=sp[i];
1525    }
1526  }
1527  if(m==0)
1528  {
1529    for(i=n;i>=1;i--)
1530    {
1531      m[i]=1;
1532    }
1533  }
1534
1535  int k;
1536  ideal a0;
1537  intvec m0;
1538  list V0;
1539  number a1;
1540  int m1;
1541  for(i=n;i>=1;i--)
1542  {
1543    if(m[i]!=0)
1544    {
1545      for(j=i-1;j>=1;j--)
1546      {
1547        if(m[j]!=0)
1548        {
1549          if(number(a[i])>number(a[j]))
1550          {
1551            a1=number(a[i]);
1552            a[i]=a[j];
1553            a[j]=a1;
1554            m1=m[i];
1555            m[i]=m[j];
1556            m[j]=m1;
1557            if(size(V)>0)
1558            {
1559              v=V[i];
1560              V[i]=V[j];
1561              V[j]=v;
1562            }
1563          }
1564          if(number(a[i])==number(a[j]))
1565          {
1566            m[i]=m[i]+m[j];
1567            m[j]=0;
1568            if(size(V)>0)
1569            {
1570              V[i]=V[i]+V[j];
1571            }
1572          }
1573        }
1574      }
1575      k++;
1576      a0[k]=a[i];
1577      m0[k]=m[i];
1578      if(size(V)>0)
1579      {
1580        V0[k]=V[i];
1581      }
1582    }
1583  }
1584
1585  if(size(V0)>0)
1586  {
1587    n=size(V0);
1588    module U=std(V0[n]);
1589    for(i=n-1;i>=1;i--)
1590    {
1591      V0[i]=simplify(reduce(V0[i],U),1);
1592      if(i>=2)
1593      {
1594        U=std(U+V0[i]);
1595      }
1596    }
1597  }
1598
1599  if(k>0)
1600  {
1601    sp=a0,m0;
1602    if(size(V0)>0)
1603    {
1604      sp[3]=V0;
1605    }
1606  }
1607  return(sp);
1608}
1609example
1610{ "EXAMPLE:"; echo=2;
1611  ring R=0,(x,y),ds;
1612  list sp=list(ideal(-1/2,-3/10,-3/10,-1/10,-1/10,0,1/10,1/10,3/10,3/10,1/2));
1613  spprint(spnf(sp));
1614}
1615///////////////////////////////////////////////////////////////////////////////
1616
1617proc spprint(list sp)
1618"USAGE:   spprint(sp); list sp (helper routine for spnf)
1619RETURN:  string s;  spectrum sp
1620EXAMPLE: example spprint; shows examples
1621SEE ALSO: gmssing_lib, spnf
1622"
1623{
1624  string s;
1625  for(int i=1;i<size(sp[2]);i++)
1626  {
1627    s=s+"("+string(sp[1][i])+","+string(sp[2][i])+"),";
1628  }
1629  s=s+"("+string(sp[1][i])+","+string(sp[2][i])+")";
1630  return(s);
1631}
1632example
1633{ "EXAMPLE:"; echo=2;
1634  ring R=0,(x,y),ds;
1635  list sp=list(ideal(-1/2,-3/10,-1/10,0,1/10,3/10,1/2),intvec(1,2,2,1,2,2,1));
1636  spprint(sp);
1637}
1638///////////////////////////////////////////////////////////////////////////////
1639
1640proc jordan(matrix M,list #)
1641"USAGE:   jordan(M); matrix M
1642ASSUME:  eigenvalues of M in basefield
1643RETURN:
1644@format
1645list l;  Jordan data of M
1646  ideal l[1];
1647    number l[1][i];  eigenvalue of i-th Jordan block of M
1648  intvec l[2];
1649    int l[2][i];  size of i-th Jordan block of M
1650  intvec l[3];
1651    int l[3][i];  multiplicity of i-th Jordan block of M
1652@end format
1653EXAMPLE: example jordan; shows examples
1654"
1655{
1656  if(nrows(M)==0)
1657  {
1658    ERROR("non empty expected");
1659  }
1660  if(ncols(M)!=nrows(M))
1661  {
1662    ERROR("square matrix expected");
1663  }
1664
1665  M=jet(M,0);
1666
1667  if(size(#)==0)
1668  {
1669    #=eigenvals(M);
1670  }
1671  def e0,m0=#[1..2];
1672
1673  int i;
1674  for(i=1;i<=ncols(e0);i++)
1675  {
1676    if(deg(e0[i])>0)
1677    {
1678
1679      ERROR("eigenvalues in coefficient field expected");
1680      return(list());
1681    }
1682  }
1683
1684  int j,k;
1685  matrix N0,N1;
1686  module K0;
1687  list K;
1688  ideal e;
1689  intvec s,m;
1690
1691  for(i=1;i<=ncols(e0);i++)
1692  {
1693    N0=M-e0[i]*matrix(freemodule(ncols(M)));
1694
1695    N1=N0;
1696    K0=0;
1697    K=module();
1698    while(size(K0)<m0[i])
1699    {
1700      K0=syz(N1);
1701      K=K+list(K0);
1702      N1=N1*N0;
1703    }
1704
1705    for(j=2;j<size(K);j++)
1706    {
1707      if(2*size(K[j])-size(K[j-1])-size(K[j+1])>0)
1708      {
1709        k++;
1710        e[k]=e0[i];
1711        s[k]=j-1;
1712        m[k]=2*size(K[j])-size(K[j-1])-size(K[j+1]);
1713      }
1714    }
1715    if(size(K[j])-size(K[j-1])>0)
1716    {
1717      k++;
1718      e[k]=e0[i];
1719      s[k]=j-1;
1720      m[k]=size(K[j])-size(K[j-1]);
1721    }
1722  }
1723
1724  return(list(e,s,m));
1725}
1726example
1727{ "EXAMPLE:"; echo=2;
1728  ring R=0,x,dp;
1729  matrix M[3][3]=3,2,1,0,2,1,0,0,3;
1730  print(M);
1731  jordan(M);
1732}
1733///////////////////////////////////////////////////////////////////////////////
1734
1735proc jordanbasis(matrix M,list #)
1736"USAGE:   jordanbasis(M); matrix M
1737ASSUME:  eigenvalues of M in basefield
1738RETURN:
1739@format
1740list l:
1741  module l[1];  inverse(l[1])*M*l[1] in Jordan normal form
1742  intvec l[2];
1743    int l[2][i];  weight filtration index of l[1][i]
1744@end format
1745EXAMPLE: example jordanbasis; shows examples
1746"
1747{
1748  if(nrows(M)==0)
1749  {
1750    ERROR("non empty matrix expected");
1751  }
1752  if(ncols(M)!=nrows(M))
1753  {
1754    ERROR("square matrix expected");
1755  }
1756
1757  M=jet(M,0);
1758
1759  if(size(#)==0)
1760  {
1761    #=eigenvals(M);
1762  }
1763  def e,m=#[1..2];
1764
1765  for(int i=1;i<=ncols(e);i++)
1766  {
1767    if(deg(e[i])>0)
1768    {
1769      ERROR("eigenvalues in coefficient field expected");
1770      return(freemodule(ncols(M)));
1771    }
1772  }
1773
1774  int j,k,l,n;
1775  matrix N0,N1;
1776  module K0,K1;
1777  list K;
1778  matrix u[ncols(M)][1];
1779  module U;
1780  intvec w;
1781
1782  for(i=1;i<=ncols(e);i++)
1783  {
1784    N0=M-e[i]*matrix(freemodule(ncols(M)));
1785
1786    N1=N0;
1787    K0=0;
1788    K=list();
1789    while(size(K0)<m[i])
1790    {
1791      K0=syz(N1);
1792      K=K+list(K0);
1793      N1=N1*N0;
1794    }
1795
1796    K1=0;
1797    for(j=1;j<size(K);j++)
1798    {
1799      K0=K[j];
1800      K[j]=interred(reduce(K[j],std(K1+module(N0*K[j+1]))));
1801      K1=K0;
1802    }
1803    K[j]=interred(reduce(K[j],std(K1)));
1804
1805    for(l=size(K);l>=1;l--)
1806    {
1807      for(k=size(K[l]);k>0;k--)
1808      {
1809        u=K[l][k];
1810        for(j=l;j>=1;j--)
1811        {
1812          U=U+module(u);
1813          n++;
1814          w[n]=2*j-l-1;
1815          u=N0*u;
1816        }
1817      }
1818    }
1819  }
1820
1821  return(list(U,w));
1822}
1823example
1824{ "EXAMPLE:"; echo=2;
1825  ring R=0,x,dp;
1826  matrix M[3][3]=3,2,1,0,2,1,0,0,3;
1827  print(M);
1828  list l=jordanbasis(M);
1829  print(l[1]);
1830  print(l[2]);
1831  print(inverse(l[1])*M*l[1]);
1832}
1833///////////////////////////////////////////////////////////////////////////////
1834
1835proc jordanmatrix(list jd)
1836"USAGE:   jordanmatrix(list(e,s,m)); ideal e, intvec s, intvec m
1837ASSUME:  ncols(e)==size(s)==size(m)
1838RETURN:
1839@format
1840matrix J;  Jordan matrix with list(e,s,m)==jordan(J)
1841@end format
1842EXAMPLE: example jordanmatrix; shows examples
1843"
1844{
1845  ideal e=jd[1];
1846  intvec s=jd[2];
1847  intvec m=jd[3];
1848  if(ncols(e)!=size(s)||ncols(e)!=size(m))
1849  {
1850    ERROR("arguments of equal size expected");
1851  }
1852
1853  int i,j,k,l;
1854  int n=int((transpose(matrix(s))*matrix(m))[1,1]);
1855  matrix J[n][n];
1856  for(k=1;k<=ncols(e);k++)
1857  {
1858    for(l=1;l<=m[k];l++)
1859    {
1860      j++;
1861      J[j,j]=e[k];
1862      for(i=s[k];i>=2;i--)
1863      {
1864        J[j+1,j]=1;
1865        j++;
1866        J[j,j]=e[k];
1867      }
1868    }
1869  }
1870
1871  return(J);
1872}
1873example
1874{ "EXAMPLE:"; echo=2;
1875  ring R=0,x,dp;
1876  ideal e=ideal(2,3);
1877  intvec s=1,2;
1878  intvec m=1,1;
1879  print(jordanmatrix(list(e,s,m)));
1880}
1881///////////////////////////////////////////////////////////////////////////////
1882
1883proc jordannf(matrix M,list #)
1884"USAGE:   jordannf(M); matrix M
1885ASSUME:  eigenvalues of M in basefield
1886RETURN:  matrix J; Jordan normal form of M
1887EXAMPLE: example jordannf; shows examples
1888"
1889{
1890  return(jordanmatrix(jordan(M,#)));
1891}
1892example
1893{ "EXAMPLE:"; echo=2;
1894  ring R=0,x,dp;
1895  matrix M[3][3]=3,2,1,0,2,1,0,0,3;
1896  print(M);
1897  print(jordannf(M));
1898}
1899
1900///////////////////////////////////////////////////////////////////////////////
1901
1902/*
1903///////////////////////////////////////////////////////////////////////////////
1904//          Auskommentierte zusaetzliche Beispiele
1905//
1906///////////////////////////////////////////////////////////////////////////////
1907// Singular for ix86-Linux version 1-3-10  (2000121517)  Dec 15 2000 17:55:12
1908// Rechnungen auf AMD700 mit 632 MB
1909
1910 LIB "linalg.lib";
1911
19121. Sparse integer Matrizen
1913--------------------------
1914ring r1=0,(x),dp;
1915system("--random", 12345678);
1916int n = 70;
1917matrix m = sparsemat(n,n,50,100);
1918option(prot,mem);
1919
1920int t=timer;
1921matrix im = inverse(m,1)[1];
1922timer-t;
1923print(im*m);
1924//list l0 = watchdog(100,"inverse("+"m"+",3)");
1925//bricht bei 100 sec ab und gibt l0[1]: string Killed zurueck
1926
1927//inverse(m,1): std       5sec 5,5 MB
1928//inverse(m,2): interred  12sec
1929//inverse(m,2): lift      nach 180 sec 13MB abgebrochen
1930//n=60: linalgorig: 3  linalg: 5
1931//n=70: linalgorig: 6,7 linalg: 11,12
1932// aber linalgorig rechnet falsch!
1933
19342. Sparse poly Matrizen
1935-----------------------
1936ring r=(0),(a,b,c),dp;
1937system("--random", 12345678);
1938int n=6;
1939matrix m = sparsematrix(n,n,2,0,50,50,9); //matrix of polys of deg <=2
1940option(prot,mem);
1941
1942int t=timer;
1943matrix im = inverse(m);
1944timer-t;
1945print(im*m);
1946//inverse(m,1): std       0sec 1MB
1947//inverse(m,2): interred  0sec 1MB
1948//inverse(m,2): lift      nach 2000 sec 33MB abgebrochen
1949
19503. Sparse Matrizen mit Parametern
1951---------------------------------
1952//liborig rechnet hier falsch!
1953ring r=(0),(a,b),dp;
1954system("--random", 12345678);
1955int n=7;
1956matrix m = sparsematrix(n,n,1,0,40,50,9);
1957ring r1 = (0,a,b),(x),dp;
1958matrix m = imap(r,m);
1959option(prot,mem);
1960
1961int t=timer;
1962matrix im = inverse(m);
1963timer-t;
1964print(im*m);
1965//inverse(m)=inverse(m,3):15 sec inverse(m,1)=1sec inverse(m,2):>120sec
1966//Bei Parametern vergeht die Zeit beim Normieren!
1967
19683. Sparse Matrizen mit Variablen und Parametern
1969-----------------------------------------------
1970ring r=(0),(a,b),dp;
1971system("--random", 12345678);
1972int n=6;
1973matrix m = sparsematrix(n,n,1,0,35,50,9);
1974ring r1 = (0,a),(b),dp;
1975matrix m = imap(r,m);
1976option(prot,mem);
1977
1978int t=timer;
1979matrix im = inverse(m,3);
1980timer-t;
1981print(im*m);
1982//n=7: inverse(m,3):lange sec inverse(m,1)=1sec inverse(m,2):1sec
1983
19844. Ueber Polynomring invertierbare Matrizen
1985-------------------------------------------
1986LIB"random.lib"; LIB"linalg.lib";
1987system("--random", 12345678);
1988int n =3;
1989ring r= 0,(x,y,z),(C,dp);
1990matrix A=triagmatrix(n,n,1,0,0,50,2);
1991intmat B=sparsetriag(n,n,20,1);
1992matrix M = A*transpose(B);
1993M=M*transpose(M);
1994M[1,1..ncols(M)]=M[1,1..n]+xyz*M[n,1..ncols(M)];
1995print(M);
1996//M hat det=1 nach Konstruktion
1997
1998int t=timer;
1999matrix iM=inverse(M);
2000timer-t;
2001print(iM*M);                      //test
2002
2003//ACHTUNG: Interred liefert i.A. keine Inverse, Gegenbeispiel z.B.
2004//mit n=3
2005//eifacheres Gegenbeispiel:
2006matrix M =
20079yz+3y+3z+2,             9y2+6y+1,
20089xyz+3xy+3xz-9z2+2x-6z-1,9xy2+6xy-9yz+x-3y-3z
2009//det M=1, inverse(M,2); ->// ** matrix is not invertible
2010//lead(M); 9xyz*gen(2) 9xy2*gen(2) nicht teilbar!
2011
20125. charpoly:
2013-----------
2014//ring rp=(0,A,B,C),(x),dp;
2015ring r=0,(A,B,C,x),dp;
2016matrix m[12][12]=
2017AC,BC,-3BC,0,-A2+B2,-3AC+1,B2,   B2,  1,   0, -C2+1,0,
20181, 1, 2C,  0,0,     B,     -A,   -4C, 2A+1,0, 0,    0,
20190, 0, 0,   1,0,     2C+1,  -4C+1,-A,  B+1, 0, B+1,  3B,
2020AB,B2,0,   1,0,     1,     0,    1,   A,   0, 1,    B+1,
20211, 0, 1,   0,0,     1,     0,    -C2, 0,   1, 0,    1,
20220, 0, 2,   1,2A,    1,     0,    0,   0,   0, 1,    1,
20230, 1, 0,   1,1,     2,     A,    3B+1,1,   B2,1,    1,
20240, 1, 0,   1,1,     1,     1,    1,   2,   0, 0,    0,
20251, 0, 1,   0,0,     0,     1,    0,   1,   1, 0,    3,
20261, 3B,B2+1,0,0,     1,     0,    1,   0,   0, 1,    0,
20270, 0, 1,   0,0,     0,     0,    1,   0,   0, 0,    0,
20280, 1, 0,   1,1,     3,     3B+1, 0,   1,   1, 1,    0;
2029option(prot,mem);
2030
2031int t=timer;
2032poly q=charpoly(m,"x");         //1sec, charpoly_B 1sec, 16MB
2033timer-t;
2034//1sec, charpoly_B 1sec, 16MB  (gleich in r und rp)
2035
2036*/
Note: See TracBrowser for help on using the repository browser.