source: git/Singular/LIB/curvepar.lib @ f5d2647

spielwiese
Last change on this file since f5d2647 was 3686937, checked in by Oleksandr Motsak <motsak@…>, 11 years ago
Added '$Id$' as a comment to all libs (LIB/*.lib)
  • Property mode set to 100644
File size: 39.3 KB
Line 
1////////////////////////////////////////////////////////////////////////////////
2version="version curvepar.lib 4.0.0.0 Jun_2013 "; // $Id$
3category="Singularity Theory";
4info="
5LIBRARY:  curvepar.lib  Resolution of space curve singularities, semi-group
6
7AUTHOR:     Gerhard Pfister    email: pfister@mathematik.uni-kl.de
8            Nil Sahin          email: e150916@metu.edu.tr
9            Maryna Viazovska   email: viazovsk@mathematik.uni-kl.de
10
11SEE ALSO: spcurve_lib
12
13PROCEDURES:
14BlowingUp(f,I,l);             BlowingUp of V(I) at the point 0;
15CurveRes(I);                  Resolution of V(I)
16CurveParam(I);                Parametrization of algebraic branches of V(I)
17WSemigroup(X,b);              Weierstrass semigroup of the curve
18primparam(x,y,c);             HN matrix of parametrization(x(t),y(t))
19MultiplicitySequence(I);      Multiplicity sequences of the branches of plane curve V(I)
20CharacteristicExponents(I);   Characteristic exponents of the branches of plane curve V(I)
21IntersectionMatrix(I);        Intersection Matrix of the branches of plane curve V(I)
22ContactMatrix(I);             Contact Matrix of the branches of plane curve V(I)
23plainInvariants(I);           Invariants of the branches of plane curve V(I)
24";
25
26LIB "sing.lib";
27LIB "primdec.lib";
28LIB "linalg.lib";
29LIB "ring.lib";
30LIB "alexpoly.lib";
31LIB "matrix.lib";
32
33//////////////////////////////////////////////////////////////
34//----------Resolution of singular curve--------------------//
35//////////////////////////////////////////////////////////////
36
37proc BlowingUp(poly f,ideal I,list l,list #)
38"USAGE:   BlowingUp(f,I,l);
39          f=poly
40          b=ideal
41          l=list
42
43ASSUME:   The basering is r=0,(x(1..n),a),dp
44          f is an irrreducible polynomial in k[a],
45          I is an ideal of a curve(if we consider a as a parameter)
46
47COMPUTE:  Blowing-up of the curve at point 0.
48
49RETURN:   list C of charts.
50          Each chart C[i] is a list of size 5 (reps. 6 in case of plane curves)
51          C[i][1] is an integer j. It shows, which standard chart do we consider.
52          C[i][2] is an irreducible poly g in k[a]. It is a minimal polynomial
53                  for the new parameter.
54          C[i][3] is an ideal H in k[a].
55                  c_i=F_i(a_new) for i=1..n,
56                  a_old=H[n+1](a_new).
57          C[i][4] is a map teta:k[x(1)..x(n),a]-->k[x(1)..x(n),a] from the new
58                  curve to the old one.
59                  x(1)-->x(j)*x(1)
60                  . . .
61                  x(j)-->x(j)
62                  . . .
63                  x(n)-->x(j)*(c_n+x(n))
64          C[i][5] is an ideal J of a new curve. J=teta(I).
65          C[i][6] is the list of exceptional divisors in the chart
66
67EXAMPLE: example BlowingUp; shows an example"
68{
69  def r=basering;
70  int n=nvars(r)-1;
71  ring r1=(0,a),(x(1..n)),ds;
72  number f=leadcoef(imap(r,f));
73  minpoly=f;
74  ideal I=imap(r,I);
75  ideal locI=std(I);
76  ideal J=tangentcone(I);
77  setring r;
78  ideal J=imap(r1,J);
79  ideal locI=imap(r1,locI);
80  int j;
81  int i;
82  list C,E;
83  list C1;
84  ideal B;
85  poly g;
86  ideal F;
87  poly b,p;
88  list Z;
89  list Z1;
90  ideal D;
91  map teta;
92  ideal D1;
93  map teta1;
94  int k,e;
95  ideal I1;
96  ideal I2;
97  int ind;
98  list w=mlist(l,n);
99  for(j=1;j<=n;j++)
100  {
101    B=J;
102    for(i=1;i<j;i++) {B=B+x(w[i]);}
103    B=B+(x(w[j])-1);
104    B=B+f;
105    Z=Z1;
106    if(dim(std(B))==0)
107    {
108      Z=ZeroIdeal(B);
109      for(i=1;i<j;i++)
110      {
111        D[w[i]]=x(w[j])*x(w[i]);
112      }
113      D[w[j]]=x(w[j]);
114      for(i=j+1;i<=n;i++)
115      {
116        D[w[i]]=x(w[j])*x(w[i]);
117      }
118      D[n+1]=a;
119      teta=r,D;
120      I1=teta(locI);
121      I1=reduce(I1,std(f));
122      ind=1;
123      for(i=1;i<=size(I1);i++)
124      {
125        ind=1;
126        while(ind==1)
127        {
128          if(gcd(x(w[j]),I1[i])==x(w[j])){I1[i]=I1[i]/x(w[j]);}
129          else{ind=0;}
130        }
131      }
132    }
133    for(k=1;k<=size(Z);k++)
134    {
135      g=Z[k][1];
136      for(i=1;i<=n;i++){F[i]=Z[k][2][i];}
137      b=Z[k][3];
138      C1[1]=w[j];
139      C1[2]=g;
140      C1[3]=F;
141      for(i=1;i<j;i++)
142      {D[w[i]]=x(w[j])*x(w[i]);}
143      D[w[j]]=x(w[j]);
144      for(i=j+1;i<=n;i++)
145      {D[w[i]]=x(w[j])*(F[w[i]]+x(w[i]));}
146      D[n+1]=Z[k][2][n+1];
147      teta=r,D;
148      C1[4]=D;
149      for(i=1;i<=j;i++)
150      {D1[w[i]]=x(w[i]);}
151      for(i=j+1;i<=n;i++)
152      {D1[w[i]]=F[w[i]]+x(w[i]);}
153      D1[n+1]=a;
154      teta1=r,D1;
155      if(nvars(basering)==3)
156      {
157        I2=quickSubst(I1[1],teta1[1],teta1[2],std(g));
158      }
159      else
160      {
161        I2=teta1(I1);
162        I2=reduce(I2,std(g));
163      }
164      C1[5]=I2;
165      if(nvars(basering)==3)
166      {
167        if(size(#)>0)
168        {
169          E=#;
170          E=teta(E);
171          for(e=1;e<=size(E);e++)
172          {
173            p=E[e];
174            while(subst(p,x(w[j]),0)==0)
175            {
176              p=p/x(w[j]);
177            }
178            if((deg(E[e])>0)&&(deg(p)==0))
179            {
180              E[e]=size(E);
181            }
182            else
183            {
184              E[e]=p;
185            }
186          }
187          E[size(E)+1]=x(w[j]);
188          C1[6]=E;
189        }
190        else
191        {
192          C1[6]=list(x(w[j]));
193        }
194      }
195      C=insert(C,C1);
196    }
197  }
198  return(C);
199}
200example
201{
202  "EXAMPLE:";echo = 2;
203  ring r=0,(x(1..3),a),dp;
204  poly f=a2+1;
205  ideal i=x(1)^2+a*x(2)^3,x(3)^2-x(2);
206  list l=1,3,2;
207  list B=BlowingUp(f,i,l);
208  B;
209}
210//============= ACHTUNG ZeroIdeal ueberarbeiten / minAssGTZ rein  ========================
211//////////////////////////////////////////////////////////////////////////////////////////
212static proc ZeroIdeal(ideal J)
213
214"USAGE:   ZeroIdeal(J);
215          J=ideal
216
217ASSUME:   J is a zero-dimensional ideal in k[x(1),...,x(n)].
218
219COMPUTE:  Primary decomposition of radical(J). Each prime ideal J[i] has the form:
220          x(1)-f[1](b),...,x(n)-f[n](b),
221          f(b)=0, f irreducible
222          for some b=x(1)*a(1)+...+x(n)*a(n), a(i) in k.
223
224RETURN:   list Z of lists.
225          Each list Z[k] is a list of size 3
226          Z[k][1] is a poly f(b)
227          Z[k][2] is an ideal H, H[n]=f[n],
228          Z[k][3] is a poly x(1)*a(1)+...+x(n)*a(n)
229
230EXAMPLE:"
231{
232  intvec opt = option(get);
233  def r=basering;
234  int n=nvars(r);
235  if(dim(std(J))!=0){return(0);}
236  ring s=0,(x(1..n)),lp;
237  ideal A; ideal S; int i; int j;
238  for(i=1;i<=n;i++) {A[i]=x(i);}
239  map phi=r,A;
240  ideal J=phi(J);
241  ideal I=radical(J);
242  list D=zerodec(I);
243  list Z; ideal H; intvec w; intvec v; int ind; ideal T; map tau; int q; list u;
244  ideal Di; poly h;
245  for(i=1;i<=size(D);i++)
246  {
247    option(redSB);
248    ind=0;q=n;
249    while(ind==0 and q>0)
250    {
251      for(j=1;j<=n;j++){T[j]=x(j);}
252      T[q]=x(n);
253      T[n]=x(q);
254      tau=s,T;
255      Di=D[i];
256      S=std(tau(Di));
257      ind=1;
258      v=leadexp(S[1]);
259      if(leadmonom(S[1])!=x(n)^v[n]){ind=0;}
260      for(j=2;j<=n;j++)
261      {
262        if(leadmonom(S[j])!=x(n-j+1)){ind=0;}
263        H[n-j+1]= -S[j]/leadcoef(S[j])+x(n-j+1);
264        v=leadexp(H[n-j+1]);
265        if(leadcoef(H[n-j+1])*leadmonom(H[n-j+1])!=leadcoef(H[n-j+1])*x(n)^v[n])
266        {ind=0;}
267      }
268      if(ind==1)
269      {
270        u[1]=S[1];
271        H[n]=x(n);
272        H[n]=H[q];
273        H[q]=x(n);
274        u[2]=H;
275        u[3]=x(q);
276        Z[i]=u;
277      }
278      q--;
279    }
280    if(ind==0)
281    {
282      vector a;
283      while(ind==0)
284      {
285        h=x(n);
286        for(j=1;j<=n-1;j++){a=a+random(-10,10)*gen(j);h=h+a[j]*x(j);}
287        T=subst(S,x(n),h);
288        option(redSB);
289        T=std(T);
290        ind=1;
291        w=leadexp(T[1]);
292        if(leadmonom(T[1])!=x(n)^w[n]){ind=0;}
293        for(j=2;j<=n;j++)
294        {
295          if(leadmonom(T[j])!=x(n-j+1)){ind=0;}
296          H[n-j+1]= -T[j]/leadcoef(T[j])+x(n-j+1);
297          w=leadexp(H[n-j+1]);
298          if(leadmonom(H[n-j+1])*leadcoef(H[n-j+1])!=leadcoef(H[n-j+1])*x(n)^w[n])
299          {ind=0;}
300        }
301        if(ind==1)
302        {
303          list l;
304          l[1]=T[1];
305          H[n]=x(n);h=x(n);
306          for(j=1;j<=n-1;j++){H[n]=H[n]+a[j]*H[j];h=h-a[j]*x(j);}
307          l[2]=H;
308          l[3]=h;
309          Z[i]=l;
310        }
311      }
312    }
313  }
314  setring r;
315  ideal A;
316  list Z;
317  for(i=1;i<=n;i++)
318  {A[i]=var(i);}
319  map psi=s,A;
320  Z=psi(Z);
321  option(set, opt);
322  return(Z);
323}
324/////////////////////////////////////////////////////////////////////////////////////////////
325//assume that the basering is k[x(1),...,x(n),a]
326
327static proc main(ideal I,ideal Psi,poly f,list m,list l,list HN,intvec v,list HI,list #)
328{
329  def s=basering;
330  int i,z;
331  int j;
332  list C,E,resTree;
333  list C1;
334  list C2;
335  list C3;
336  list l1;
337  C2[8]=HI;
338  list m1;
339  list HN1;
340  ideal J;
341  map psi;
342  intvec w;
343  z=(SmoothTest(I,f)==1);
344  if((nvars(basering)==3)&&z&&(size(#)>0))
345  {
346    z=transversalTest(I,f,#);
347  }
348  if(z)
349  {
350    C2[1]=I;
351    C2[2]=Psi;
352    C2[3]=f;
353    C2[4]=m;
354    C2[5]=l;
355    C2[6]=HN;
356    if(nvars(basering)==3)
357    {
358      if(size(#)>0)
359      {
360        C2[9]=#;
361      }
362      C2[7]=v;
363    }
364    //C2[8][size(C2[8])+1]=list(C2[7],C2[9]);
365    C[1]=C2;
366  }
367  if(!z)
368  {
369    int mm=mmult(I,f);
370    m1=insert(m,mm,size(m));
371    if(nvars(basering)==3)
372    {
373      if(size(#)>0)
374      {
375        E=#;
376        C1=BlowingUp(f,I,l,E);
377      }
378      else
379      {
380        C1=BlowingUp(f,I,l);
381      }
382    }
383    else
384    {
385      C1=BlowingUp(f,I,l);
386    }
387    for(j=1;j<=size(C1);j++)
388    {
389      C2[1]=C1[j][5];
390      J=C1[j][4];
391      psi=s,J;
392      C2[2]=psi(Psi);
393      C2[3]=C1[j][2];
394      C2[4]=m1;
395      l1=insert(l,C1[j][1],size(l));
396      C2[5]=l1;
397      HN1=psi(HN);
398      HN1=insert(HN1,C1[j][3],size(HN1)-1);
399      C2[6]=HN1;
400      if(deg(C2[3])>1)
401      {
402        w=v,-j;
403      }
404      else
405      {
406        w=v,j;
407      }
408      C2[7]=w;
409      if(nvars(basering)==3)
410      {
411        C2[9]=C1[j][6];
412        C2[8][size(C2[8])+1]=list(C2[7],C2[9]);
413        C3=main(C2[1],C2[2],C2[3],C2[4],C2[5],C2[6],C2[7],C2[8],C2[9]);
414        C=C+C3;
415      }
416      else
417      {
418        C3=main(C2[1],C2[2],C2[3],C2[4],C2[5],C2[6],C2[7],C2[8]);
419        C=C+C3;
420      }
421    }
422  }
423  return(C);
424}
425////////////////////////////////////////////////////////////////////////////////////////////////
426
427static proc transversalTest(ideal I,poly f,list L)
428{
429  def r=basering;
430  int n=nvars(r)-1;
431  int i;
432  ring r1=(0,a),(x(1..n)),ds;
433  number f=leadcoef(imap(r,f));
434  minpoly=f;
435  ideal I=imap(r,I);
436  list L=imap(r,L);
437  ideal K=jet(L[size(L)],deg(lead(L[size(L)])));
438  ideal T=1;
439  if(size(L)>1)
440  {
441    for(i=1;i<=size(L)-1;i++)
442    {
443      if(subst(L[i],x(1),0,x(2),0)==0) break;
444    }
445    if(i<=size(L)-1)
446    {
447      T=jet(L[i],deg(lead(L[i])));
448    }
449  }
450  ideal J=jet(I[1],deg(lead(I[1])));
451  setring r;
452  ideal J=imap(r1,J);
453  ideal K=imap(r1,K);
454  ideal T=imap(r1,T);
455  int m=size(reduce(J,std(K)))+size(reduce(K,std(J)));
456  if(m)
457  {
458    m=size(reduce(J+K+T,std(ideal(x(1),x(2)))));
459  }
460  return(m);
461}
462////////////////////////////////////////////////////////////////////////////////////////////////
463static proc SmoothTest(ideal I,poly f)
464//Assume I is a radical ideal of dimension 1 in a ring k[x(1..n),a]
465//Returns 1 if a curve V(I) is smooth at point 0 and returns 0 otherwise
466{
467  int ind;
468  int l;
469  def t=basering;
470  int n=nvars(t)-1;
471  ring r1=(0,a),(x(1..n)),dp;
472  number f=leadcoef(imap(t,f));
473  minpoly=f;
474  ideal I=imap(t,I);
475  matrix M=jacob(I);
476  for(l=1;l<=n;l++){M=subst(M,x(l),0);}
477  if(mat_rk(M)==(n-1)){ind=1;}
478  return(ind);
479}
480////////////////////////////////////////////////////////////////////////////////////////////////
481proc CurveRes(ideal I)
482"USAGE:   CurveRes(I);
483          I ideal
484ASSUME:   The basering is r=0,(x(1..n))
485          V(I) is a curve with a singular point 0.
486COMPUTE:  Resolution of the curve V(I).
487RETURN:   a ring R=basering+k[a]
488          Ring R contains a list Resolve
489          Resolve is a list of charts
490          Each Resolve[i] is a list of size 6
491          Resolve[i][1] is an ideal J of a new curve. J=teta(I).
492          Resolve[i][2] ideal which represents the map
493                        teta:k[x(1)..x(n),a]-->k[x(1)..x(n),a] from the
494                        new curve to the old one.
495          Resolve[i][3] is an irreducible poly g in k[a]. It is a minimal polynomial for the
496                        new parameter a. deg(g) gives the number of branches in Resolve[i]
497          Resolve[i][4] sequence of multiplicities (sum over all branches in Resolve as long as
498                        they intersect each other !)
499          Resolve[i][5] is a list of integers l. It shows, which standard charts we considered.
500          Resolve[i][6] HN matrix
501          Resolve[i][7] (only for plane curves) the development of exceptional divisors
502                        the entries correspond to the i-th blowing up. The first entry is an
503                        intvec. The first negative entry gives the splitting of the (over Q
504                        irreducible) branches. The second entry is a list of the exceptional
505                        divisors. If the entry is an integer i, it says that the divisor is not
506                        visible in this chart after the i-th blowing up.
507
508EXAMPLE: example CurveRes; shows an example"
509{
510  def r=basering;
511  int n=nvars(r);
512  ring s=0,(x(1..n),a),dp;
513  ideal A;
514  int i;
515  int j;
516  for(i=1;i<=n;i++){A[i]=x(i);}
517  map phi=r,A;
518  ideal I=phi(I);
519  poly f=a;
520  list l;
521  list m;
522  list HN=x(1);
523  ideal psi;
524  for(i=1;i<=n;i++){psi[i]=x(i);}
525  psi[n+1]=a;
526  intvec v;
527  list L,Resolve;
528  if(n==2)
529  {
530    ideal J=factorize(I[1],1);
531    list resolve;
532    for(int k=1;k<=size(J);k++)
533    {
534      I=J[k];
535      resolve=main(I,psi,f,m,l,HN,v,L);
536      for(i=1;i<=size(resolve);i++)
537      {
538        resolve[i][6]=delete(resolve[i][6],size(resolve[i][6]));
539        if(size(resolve[i])>=9){resolve[i]=delete(resolve[i],9);}
540        resolve[i]=delete(resolve[i],7);
541      }
542      if(k==1){Resolve=resolve;}
543      else{Resolve=Resolve+resolve;}
544    }
545  }
546  else
547  {
548    Resolve=main(I,psi,f,m,l,HN,v,L);
549    for(i=1;i<=size(Resolve);i++)
550    {
551      Resolve[i][6]=delete(Resolve[i][6],size(Resolve[i][6]));
552      Resolve[i]=delete(Resolve[i],8);
553    }
554  }
555  export(Resolve);
556  return(s);
557}
558example
559{
560  "EXAMPLE:"; echo=2;
561  ring r=0,(x,y,z),dp;
562  ideal i=x2-y3,z2-y5;
563  def s=CurveRes(i);
564  setring s;
565  Resolve;
566}
567//////////////////////////////////////////////////////////////////
568static proc mlist(list l,int n)
569{
570  list N;
571  list M;
572  int i;
573  int j;
574  for(i=1;i<=n;i++) {M[i]=i;}
575  N=l+M;
576  for(i=1;i<=size(N)-1;i++)
577  {
578    j=i+1;
579    while(j<=size(N))
580    {
581      if(N[i]==N[j]){N=delete(N,j);}
582      else
583      {j++;}
584    }
585  }
586  return(N);
587}
588/////////////////////////////////////////////////////////////////////
589//Assume that the basering is k[x(1..n),a]
590
591static proc mmult(ideal I,poly f)
592{
593  def r=basering;
594  int n=nvars(r)-1;
595  ring r1=(0,a),(x(1..n)),ds;
596  number f=leadcoef(imap(r,f));
597  minpoly=f;
598  ideal I=imap(r,I);
599  int m=mult(std(I));
600  return(m);
601}
602//////////////////////////////////////////////////////////////
603//--------Parametrization of smooth curve-------------------//
604//////////////////////////////////////////////////////////////
605
606////////////////////////////////////////////////////////////////////////
607//computes jacobian matrix, considering x(1..n) as variables and a(1..m) as parameters
608
609static proc mjacob(ideal I)
610{
611  def r=basering;
612  int n=nvars(r);
613  int k=size(I);
614  matrix M[k][n];
615  int i;
616  int j;
617  int l;
618  for(i=1;i<=k;i++)
619  {
620    for(j=1;j<=n;j++)
621    {
622      M[i,j]=diff(I[i],x(j));
623      for(l=1;l<=n;l++){M[i,j]=subst(M[i,j],x(l),0);}
624    }
625  }
626  return(M);
627}
628//////////////////////////////////////////////////////////
629static proc mmi(matrix M,int n)
630{
631  ideal l;
632  int k=nrows(M);
633  int i;
634  int j;
635  for(i=1;i<=k;i++)
636  {
637    l[i]=0;
638    for(j=1;j<=n;j++)
639    {
640      l[i]=l[i]+x(j)*M[i,j];
641    }
642  }
643  l=std(l);
644  int t=size(l);
645  i=1;
646  int mi=0;
647  while( mi==0 and i<=n-1)
648  {
649    if(diff(l[i],x(n-i))!=0){mi=n-i+1;}
650    else{i++;}
651  }
652  if(mi==0){mi=1;}
653  matrix Mi[k][n-1];
654  for(i=1;i<=k;i++)
655  {
656    for(j=1;j<=mi-1;j++)
657    {
658      Mi[i,j]=M[i,j];
659    }
660    for(j=mi;j<=n-1;j++)
661    {
662      Mi[i,j]=M[i,j+1];
663    }
664  }
665  list lmi=mi,Mi;
666  return(lmi);
667}
668//////////////////////////////////////////////////////////
669static proc mC(matrix Mi,int n)
670{
671  int k=nrows(Mi);
672  ideal c;
673  int i,j;
674  for(i=1;i<=n-1;i++)
675  {
676    c[i]=0;
677    for(j=1;j<=k;j++)
678    {
679      c[i]=c[i]+y(j)*Mi[j,i];
680    }
681  }
682  c=std(c);
683  return(c);
684}
685//////////////////////////////////////////////////////////
686static proc mmF(ideal C, matrix Mi,int n,int k)
687{
688  int s=size(C);
689  intvec mf;
690  int p=0;
691  int t=0;
692  int i;
693  int j;
694  int v=0;
695  for(i=s;i>0;i--)
696  {
697    p=t;
698    j=1;
699    while(t==p and p+j<=k)
700    {
701      if(diff(C[i],y(p+j))==0){j++;}
702      if(diff(C[i],y(p+j))!=0){t=p+j;v++;mf[v]=t;}
703    }
704  }
705  matrix B[n-1][n-1];
706  for(i=1;i<=n-1;i++)
707  {
708    for(j=1;j<=n-1;j++)
709    {
710      B[i,j]=Mi[ mf[i],j];
711    }
712  }
713  list mmf=mf,B;
714  return(mmf);
715}
716/////////////////////////////////////////////////////
717static proc cparam(ideal I,poly f,int n,int m,int N)
718{
719  def r=basering;
720  ring s=(0,a),(x(1..n)),lp;
721  number f=leadcoef(imap(r,f));
722  minpoly=f;
723  ideal I=imap(r,I);
724  matrix M=mjacob(I);
725  list l0=mmi(M,n);
726  int mi=l0[1];
727  matrix Mi=l0[2];
728  int k=nrows(Mi);
729  ring q=(0,a),(y(1..k)),lp;
730  number f=leadcoef(imap(r,f));
731  minpoly=f;
732  matrix Mi=imap(s,Mi);
733  ideal D=mC(Mi,n);
734  list l1=mmF(D,Mi,n,k);
735  intvec mf=l1[1];
736  matrix B=l1[2];
737  setring s;
738  matrix B=imap(q,B);
739  matrix C=inverse(B);
740  int i;
741  int j;
742  ideal P;
743  for(i=1;i<mi;i++){P[i]=x(i);}
744  P[mi]=x(n);
745  for(i=1;i<=n-mi;i++){P[mi+i]=x(mi+i-1);}
746  map phi=s,P;
747  ideal I1=phi(I);
748  if(nvars(basering)==2)
749  {
750    setring r;
751    ideal I1=imap(s,I1);
752    matrix C=imap(s,C);
753    list X;
754    matrix d[n-1][1];
755    matrix b[n-1][1];
756    ideal Q;
757    map psi;
758    int l;
759    for(i=1;i<=N;i++)
760    {
761      for(j=1;j<=n-1;j++)
762      {
763        d[j,1]=diff(I1[mf[j]],x(n));
764        for(l=1;l<=n;l++)
765        {
766          d[j,1]=subst(d[j,1],x(l),0);
767        }
768      }
769      b=-C*d;
770      I1=jet(I1,N-i+2);
771      X[i]=b;
772      for(j=1;j<=n-1;j++){Q[j]=x(n)*(b[j,1]+x(j));}
773      Q[n]=x(n);
774      I1[1]=quickSubst(I1[1],Q[1],Q[2],std(f));
775      I1=I1/x(n);
776    }
777    list Y=X,mi;
778    return(Y);
779  }
780  list X;
781  matrix d[n-1][1];
782  matrix b[n-1][1];
783  ideal Q;
784  map psi;
785  int l;
786  for(i=1;i<=N;i++)
787  {
788    for(j=1;j<=n-1;j++)
789    {
790      d[j,1]=diff(I1[mf[j]],x(n));
791      for(l=1;l<=n;l++){d[j,1]=subst(d[j,1],x(l),0);}
792    }
793    b=-C*d;
794    I1=jet(I1,N-i+2);
795    X[i]=b;
796    for(j=1;j<=n-1;j++){Q[j]=x(n)*(b[j,1]+x(j));}
797    Q[n]=x(n);
798    psi=s,Q;
799    I1=psi(I1);
800    I1=I1/x(n);
801  }
802  list Y=X,mi,var(1);
803  setring r;
804  list Y=imap(s,Y);
805  Y=delete(Y,3);
806  return(Y);
807}
808//////////////////////////////////////////////////////////////
809//--------Parametrization of singular curve-----------------//
810//////////////////////////////////////////////////////////////
811proc CurveParam (list #)
812"USAGE:   CurveParam(I);
813          I ideal
814ASSUME:   I is an ideal of a curve C with a singular point 0.
815COMPUTE:  Parametrization for algebraic branches of the curve C.
816RETURN:   list L of size 1.
817          L[1] is a ring ring rt=0,(t,a),ds;
818          Ring R contains a list Param
819          Param is a list of algebraic branches
820          Each Param[i] is a list of size 3
821          Param[i][1] is a list of polynomials
822          Param[i][2] is an irredusible polynomial f\in k[a].It is a minimal polynomial for
823                      the parameter a.
824          Param[i][3] is an integer b--upper bound for the conductor of Weierstrass semigroup
825
826EXAMPLE: example curveparam; shows an example"
827{
828  int i;
829  int j;
830  if(typeof(#[1])=="ideal")
831  {
832    int d=deg(#[1][1]);
833    def s=CurveRes(#[1]);
834  }
835  else
836  {
837    def s=#[1];
838  }
839  setring s;
840  int n=nvars(s)-1;
841  list Param;
842  list l;
843  ideal D,P,Q,T;
844  poly f;
845  map tau;
846  list Z;
847  list Y;
848  list X;
849  int mi;
850  int b;
851  int k;
852  int dd;
853  for(j=1;j<=size(Resolve);j++)
854  {
855    b=0;
856    for(k=1;k<=size(Resolve[j][4]);k++)
857    {
858      b=b+Resolve[j][4][k]*(Resolve[j][4][k]+1);
859    }
860    if((n==2)&&(size(Resolve[j][4])==0))
861    {b=d;}
862    Y=cparam(Resolve[j][1],Resolve[j][3],n,1,b);
863    X=Y[1];
864    mi=Y[2];
865    f=Resolve[j][3];
866    for(i=1;i<mi;i++)
867    {
868      P[i]=0;
869      for(k=1;k<=b;k++){P[i]=P[i]+X[k][i,1]*(x(1)^k);}
870    }
871    P[mi]=x(1);
872    for(i=mi+1;i<=n;i++)
873    {
874      P[i]=0;
875      for(k=1;k<=b;k++){P[i]=P[i]+X[k][i-1,1]*(x(1)^k);}
876    }
877    P[n+1]=a;
878    tau=s,P;
879    T=Resolve[j][2];
880//HERE A TEST FOR dd
881    if(nvars(basering)==3)
882    {
883      dd=boundparam(P[2]);
884      if(dd==1){dd=boundparam(P[1]);}
885      P[1]=jet(P[1],dd);
886      P[2]=jet(P[2],dd);
887      Q[1]=quickSubst(T[1],P[1],P[2],std(f));
888      Q[2]=quickSubst(T[2],P[1],P[2],std(f));
889      Q[3]=a;
890    }
891    else
892    {
893      Q=tau(T);
894    }
895    for(i=1;i<=n;i++){Z[i]=jet(reduce(Q[i],std(f)),b+1);}
896    l[1]=Z;
897    l[2]=f;
898    l[3]=b;
899    Param[j]=l;
900  }
901  ring rt=0,(t,a),ds;
902  ideal D;
903  D[1]=t;
904  D[n+1]=a;
905  map delta=s,D;
906  list Param=delta(Param);
907  export(Param);
908  return(rt);
909}
910example
911{
912  "EXAMPLE:";echo=2;
913  ring r=0,(x,y,z),dp;
914  ideal i=x2-y3,z2-y5;
915  def s=CurveParam(i);
916  setring s;
917  Param;
918}
919///////////////////////////////////////////////////////////////////////////////////////////
920//----------Computation of Weierstrass Semigroup from parametrization--------------------//
921///////////////////////////////////////////////////////////////////////////////////////////
922static proc Semi(intvec G,int b)
923"USAGE:   Semi(G,b);
924          G=intvec
925          b=int
926ASSUME:   G[1]<=G[2]<=...<=G[k],
927COMPUTE:  elements of semigroup S generated by the enteries of G till the bound b.
928          For each element i of S computes the list of integer vectors v of dimension
929          k=size(G), such that g[1]*v[1]+g[2]*v[2]+...+g[k]*v[k]=i. If there exists
930          conductor  of semigroup S  c<b-n, where n is minimal element of G, then
931          computes also c+n.
932RETURN:   list M of size 2.
933          L=M[1] is a list  of size min(b,c+n).
934          L[i] is a list of integer vectors.
935          If i is not in a semigroup S than L[i] is empty.
936          M[2] is an integer =min(b,c+n)
937          M[3] minimal generators of S
938EXAMPLE:"
939{
940  list L;
941  list l;
942  int i;
943  for(i=1;i<=b;i++){L[i]=l;}
944  int k=size(G);
945  int n=G[1];
946  int j;
947  int t;
948  int q;
949  int c=0;
950  intvec w;
951  intvec v;
952  for(i=1;i<=k;i++)
953  {
954    for(j=1;j<=k;j++)
955    {
956      if(j==i){w[j]=1;}
957      else{w[j]=0;}
958    }
959    L[G[i]]=insert(L[G[i]],w);
960  }
961  list L1=L;
962  int s=0;
963  int s1=0;
964  i=1;
965  int p;
966  while(i<=b and s<n)
967  {
968    for(j=1;j<=k;j++)
969    {
970      if(i-G[j]>0)
971      {
972        if(size(L[i-G[j]])>0)
973        {
974          for(t=1;t<=size(L[i-G[j]]);t++)
975          {
976            v=L[i-G[j]][t];
977            p=1;
978            for(q=1;q<j;q++)
979            {
980              if(v[q]>0){p=0;}
981            }
982            if(p==1){v[j]=v[j]+1;L[i]=insert(L[i],v);}
983          }
984        }
985      }
986    }
987    if(size(L[i])>0){s1=1;}
988    s=s1*(s+1);
989    s1=0;
990    i++;
991  }
992  intvec Gmin;
993  int jmin=1;
994  for(j=1;j<=k;j++)
995  {
996    if(size(L[G[j]])==size(L1[G[j]]) && G[j]<i)
997    {
998      Gmin[jmin]=G[j];
999      L1[G[j]]=insert(L1[G[j]],0);
1000      jmin++;
1001    }
1002  }
1003  list M=L,i-1,Gmin;
1004  return(M);
1005}
1006///////////////////////////////////////////////////////////////////////////////////////////
1007static proc AddElem(list L,int b,int k,int g,int n)
1008"ASSUME:  L list of size b. L[i] list of integer vectors of dimension k.
1009          b=int
1010          k=int as above
1011          g=int new generator
1012          n=int. minimal generator
1013RETURN: list M
1014        M[1]=new L;
1015        M[2]=new b;"
1016{
1017  int i,j;
1018  intvec v;
1019  for(i=1;i<=b;i++)
1020  {
1021    if(size(L[i])>0)
1022    {
1023      for(j=1;j<=size(L[i]);j++)
1024      {
1025        v=L[i][j];
1026        v[k+1]=0;
1027        L[i][j]=v;
1028      }
1029    }
1030  }
1031  intvec w;
1032  w[k+1]=1;
1033  L[g]=insert(L[g],w);
1034  int s=0;
1035  int s1=0;
1036  i=1;
1037  while(i<=b and s<n)
1038  {
1039    if(i-g>0)
1040    {
1041      if(size(L[i-g])>0)
1042      {
1043        for(j=1;j<=size(L[i-g]);j++)
1044        {
1045          v=L[i-g][j];
1046          v[k+1]=v[k+1]+1;
1047          L[i]=insert(L[i],v);
1048        }
1049      }
1050    }
1051    if(size(L[i])>0){s1=1;}
1052    s=s1*(s+1);
1053    s1=0;
1054    i++;
1055  }
1056  int b1=i-1;
1057  list M=L,b1;
1058  return(M);
1059}
1060///////////////////////////////////////////////////////////////////////////////////////////
1061proc WSemigroup(list X,int b0)
1062"USAGE:    WSemigroup(X,b0);
1063           X a list of polinomials in one vaiable, say t.
1064           b0 an integer
1065COMPUTE:   Weierstrass semigroup of space curve C,which is given by a parametrization
1066           X[1](t),...,X[k](t), till the bound b0.
1067
1068ASSUME:    b0 is greater then conductor
1069RETURN:    list M of size 5.
1070           M[1]= list of integers, which are minimal generators set of the Weierstrass semigroup.
1071           M[2]=integer, conductor of the Weierstrass semigroup.
1072           M[3]=intvec, all elements of the Weierstrass semigroup till some bound b,
1073           which is greather than conductor.
1074WARNING:   works only over the ring with one variable with ordering ds
1075EXAMPLE: example WSemigroup; shows an example"
1076
1077{
1078  int k=size(X);
1079  intvec G;
1080  int i,i2;
1081  poly t=var(1);
1082  poly h;
1083  int g;
1084  for(i=1;i<=k;i++)
1085  {G[i]=ord(X[i]);}
1086  for(i=1;i<k;i++)
1087  {
1088    for(i2=i;i2<=k;i2++)
1089    {
1090      if(G[i]>G[i2])
1091      {
1092        g=G[i];G[i]=G[i2];G[i2]=g;
1093        h=X[i];X[i]=X[i2];X[i2]=h;
1094      }
1095    }
1096  }
1097  list U=Semi(G,b0);
1098  list L=U[1];
1099  int b=U[2];
1100  G=U[3];
1101  int k1=size(G);
1102  list N;
1103  list l;
1104  for(i=1;i<=b;i++){N[i]=l;}
1105  int j;
1106  for(j=b0;j>b;j--){L=delete(L,j);}
1107  poly p;
1108  int s;
1109  int e;
1110  for(i=1;i<=b;i++)
1111  {
1112    for(j=1;j<=size(L[i]);j++)
1113    {
1114      p=1;
1115      for(s=1;s<=k;s++)
1116      {
1117        for(e=1;e<=L[i][j][s];e++)
1118        {
1119          p=p*X[s];
1120          p=jet(p,b);
1121        }
1122      }
1123      N[i]=insert(N[i],p);
1124    }
1125  }
1126  int j1;
1127  int j2;
1128  list M;
1129  poly c1;
1130  poly c2;
1131  poly f;
1132  int m;
1133  int b1;
1134  ideal I;
1135  matrix C;
1136  matrix C1;
1137  int q;
1138  int i1;
1139  i=1;
1140  while(i<=b)
1141  {
1142    for(j1=2;j1<=size(N[i]);j1++)
1143    {
1144      for(j2=1;j2<j1;j2++)
1145      {
1146        c1=coeffs(N[i][j1],t)[i+1,1];
1147        c2=coeffs(N[i][j2],t)[i+1,1];
1148        f=c2*N[i][j1]-c1*N[i][j2];
1149        m=ord(f);
1150        if(m>=0)
1151        {
1152          if(size(N[m])==0)
1153          {
1154            N[m]=insert(N[m],f);
1155            if(size(L[m])==0)
1156            {
1157              M=AddElem(L,b,k,m,G[1]);
1158              L=M[1];
1159              b1=M[2];
1160              G[k1+1]=m;
1161              X[k+1]=f;
1162              N[m]=insert(N[m],f);
1163              k=k+1;
1164              k1=k1+1;
1165              if(b1<b)
1166              {
1167                for(i1=1;i1<=b1;i1++)
1168                {
1169                  for(s=1;s<=size(N[i1]);s++){N[i1][s]=jet(N[i1][s],b1);}
1170                }
1171                for(s=size(N);s>b1;s--){N=delete(N,s);}
1172                for(s=size(L);s>b1;s--){L=delete(L,s);}
1173              }
1174              b=b1;
1175            }
1176          }
1177          else
1178          {
1179            for(q=1;q<=size(N[m]);q++){I[q]=N[m][q];}
1180            I[size(N[m])+1]=f;
1181            C=coeffs(I,t);
1182            C1=gauss_col(C);
1183            if(C1[size(N[m])+1]!=0){N[m]=insert(N[m],f);}
1184          }
1185        }
1186      }
1187    }
1188    i++;
1189  }
1190  intvec S;
1191  j=1;
1192  for(i=1;i<=b;i++)
1193  {
1194    if(size(L[i])>0){S[j]=i;j++;}
1195  }
1196  U=Semi(G,b);
1197  G=U[3];
1198  list Q=G,b-G[1]+1,S;
1199  return(Q);
1200}
1201example
1202{
1203  "EXAMPLE:";echo=2;
1204  ring r=0,(t),ds;
1205  list X=t4,t5+t11,t9+2*t7;
1206  list L=WSemigroup(X,30);
1207  L;
1208}
1209////////////////////////////////////////////////////////////////////////////////////////////
1210
1211static proc quickSubst(poly h, poly r, poly s,ideal I)
1212{
1213//=== computes h(r,s) mod I for h in Q[x(1),x(2),a]
1214  attrib(I,"isSB",1);
1215  if((r==x(1))&&(s==x(2))){return(reduce(h,I));}
1216  poly q1 = 1;
1217  poly q2 = 1;
1218  poly q3 = 1;
1219  int i,j,e1,e2,e3;
1220  list L,L1,L2,L3;
1221  if(r==x(1))
1222  {
1223    matrix M=coeffs(h,x(2));
1224    L[1]=1;
1225    for(i=2;i<=nrows(M);i++)
1226    {
1227      q2 = reduce(q2*s,I);
1228      L[i]=q2;
1229    }
1230    i=1;
1231    h=0;
1232    while(i <= nrows(M))
1233    {
1234      if(M[i,1]!=0)
1235      {
1236         h=h+M[i,1]*L[i];
1237      }
1238      i++;
1239    }
1240    h=reduce(h,I);
1241    return(h);
1242  }
1243  if(s==x(2))
1244  {
1245    matrix M=coeffs(h,x(1));
1246    L[1]=1;
1247    for(i=2;i<=nrows(M);i++)
1248    {
1249      q1 = reduce(q1*r,I);
1250      L[i]=q1;
1251    }
1252    i=1;
1253    h=0;
1254    while(i <= nrows(M))
1255    {
1256      if(M[i,1]!=0)
1257      {
1258        h=h+M[i,1]*L[i];
1259      }
1260      i++;
1261    }
1262    h=reduce(h,I);
1263    return(h);
1264  }
1265  for(i=1;i<=size(h);i++)
1266  {
1267    if(leadexp(h[i])[1]>e1){e1=leadexp(h[i])[1];}
1268    if(leadexp(h[i])[2]>e2){e2=leadexp(h[i])[2];}
1269    if(leadexp(h[i])[3]>e3){e3=leadexp(h[i])[3];}
1270  }
1271  for(i = 1; i <= size(h); i++)
1272  {
1273    L[i] = list(leadcoef(h[i]),leadexp(h[i]));
1274  }
1275  L1[1]=1;
1276  L2[1]=1;
1277  L3[1]=1;
1278  for(i=1;i<=e1;i++)
1279  {
1280    q1 = reduce(q1*r,I);
1281    L1[i+1]=q1;
1282  }
1283  for(i=1;i<=e2;i++)
1284  {
1285    q2 = reduce(q2*s,I);
1286    L2[i+1]=q2;
1287  }
1288  for(i=1;i<=e3;i++)
1289  {
1290    q3 = reduce(q3*var(3),I);
1291    L3[i+1]=q3;
1292  }
1293  int m=size(L);
1294  i = 1;
1295  h = 0;
1296  while(i <= m)
1297  {
1298    h=h+L[i][1]*L1[L[i][2][1]+1]*L2[L[i][2][2]+1]*L3[L[i][2][3]+1];
1299    i++;
1300  }
1301  h=reduce(h,I);
1302  return(h);
1303}
1304
1305static proc semi2char(intvec v)
1306{
1307  intvec k=v[1..2];
1308  intvec w=v[1];
1309  int i,j,p,q;
1310  for(i=2;i<size(v);i++)
1311  {
1312    w[i]=gcd(w[i-1],v[i]);
1313  }
1314  for(i=3;i<=size(v);i++)
1315  {
1316    k[i]=v[i];
1317    for(j=2;j<i;j++)
1318    {
1319      k[i]=k[i]-(w[j-1] div w[j]-1)*v[j];
1320    }
1321  }
1322  return(k);
1323}
1324/////////////////////////////////////////////////////////////////////////////////////////////////
1325proc primparam(poly x,poly y,int c)
1326"USAGE:  MultiplicitySequence(x,y,c);  x poly, y poly, c integer
1327ASSUME:  x and y are polynomials in k(a)[t] such that (x,y) is a primitive parametrization of
1328         a plane curve branch and ord(x)<ord(y).
1329RETURN:  Hamburger-Noether Matrix of the curve branch given parametrically by (x,y).
1330EXAMPLE: example primparam;  shows an example
1331"
1332{
1333  int i,h;
1334  poly F,z;
1335  list L;
1336  while(ord(x)>1)
1337  {
1338    list v;
1339    while(ord(y)>=ord(x))
1340    {
1341      F=divide(y,x,c);
1342      if(ord(F)==0)
1343      {
1344        v=insert(v,subst(F,t,0),size(v));
1345        y=F-subst(F,t,0);
1346      }
1347      else
1348      {
1349        v=insert(v,0,size(v));
1350        y=F;
1351      }
1352    }
1353    v=insert(v,t,size(v));
1354    L=insert(L,transform(v),size(L));
1355    z=x;
1356    x=y;
1357    y=z;
1358    kill v;
1359  }
1360  if(ord(x)==1)
1361  {
1362    list v;
1363    while(i<c)
1364    {
1365      F=divide(y,x,c);
1366      if(ord(F)==0)
1367      {
1368        v=insert(v,subst(F,t,0),size(v));
1369        y=F-subst(F,t,0);
1370      }
1371      else
1372      {
1373        v=insert(v,0,size(v));
1374        y=F;
1375      }
1376      if(y==0)
1377      {
1378        v=insert(v,t,size(v));
1379        break;
1380      }
1381      i++;
1382    }
1383    L=insert(L,transform(v),size(L));
1384  }
1385  return(compose(L));
1386}
1387example
1388{
1389  "EXAMPLE:"; echo=2;
1390  ring r=(0,a),t,ds;
1391  poly x=t6;
1392  poly y=t8+t11;
1393  int c=15;
1394  primparam(x,y,c);
1395}
1396
1397//////////////////////////////////////
1398//L is a list of polynomials
1399//////////////////////////////////////
1400static proc transform(list L)
1401{
1402  matrix m2;
1403  matrix m1=matrix(L[1]);
1404  for(int i=2;i<=size(L);i++)
1405  {
1406    m2=matrix(L[i]);
1407    m1=concat(m1,m2);
1408  }
1409  return(m1);
1410}
1411/////////////////////////////////////////////////////////////////
1412//L is a list of matrices
1413///////////////////////////////////////
1414static proc compose(list L)
1415{
1416  matrix M[ncols(L[1])][1]=transpose(L[1]);
1417  for(int i=2;i<=size(L);i++)
1418  {
1419    M=concat(M,transpose(L[i]));
1420  }
1421  return(transpose(M));
1422}
1423//////////////////////////////////////////////////////////
1424static proc rduce(poly p)
1425{
1426  int n=ord(p);
1427  poly q=p/(t^n);
1428  return(q);
1429}
1430////////////////////////////////////////////////////
1431static proc divide(poly p,poly q,int c)i
1432{
1433  int n=ord(p);
1434  int m=ord(q);
1435  poly p'=rduce(p);
1436  poly q'=rduce(q);
1437  poly r=t^(n-m)*p'*jet(1,q',c);
1438  return(jet(r,c));
1439}
1440/////////////////////////////////////////////////////
1441static proc contact(list L)
1442{
1443  def M=L[1];
1444  intvec v=L[2];
1445  int s,j,i;
1446  for(i=1;i<=size(v);i++)
1447  {
1448    if(v[i]<0){v[i]=-1-v[i];}
1449    for(j=1;j<=v[i];j++)
1450    {
1451      s=s+1;
1452      if(find(string(M[i,j]),"a")!=0){return(s);}
1453    }
1454  }
1455}
1456///////////////////////////////////////////////////////////
1457static proc converter(list L)
1458{
1459def s=basering;
1460list D;
1461int i,c;
1462for(i=1;i<=size(L);i++)
1463   {D=insert(D,deg(L[i][2]),size(D));}
1464ring r=(0,a),(t),ds;
1465list L=imap(s,L);
1466poly x,y,z;
1467list A,B;
1468for(i=1;i<=size(L);i++)
1469   {A[5]=D[i];
1470   x=L[i][1][1];
1471   y=L[i][1][2];
1472   if(ord(x)<=ord(y)){A[3]=0;}
1473   else{A[3]=1;
1474       z=x;
1475       x=y;
1476       y=z;
1477       }
1478   c=bound(x,y);
1479   if(c==-1){ERROR("Bound is not enough");}
1480   A[1]=primparam(x,y,c);
1481   A[2]=lengths(A[1]);
1482   A[4]=0;
1483   B=insert(B,A,size(B));
1484   A=list();
1485   }
1486ring r1=(0,a),(x,y),ds;
1487list hne=fetch(r,B);
1488export(hne);
1489return(r1);
1490}
1491//////////////////////////////////////////////////////////
1492static proc intermat(list L)
1493{
1494    int s=size(L);
1495    intvec v=L[1][5];
1496    intvec w1;
1497    int i,j,d,b,l,k,c,o,p;
1498    for(i=2;i<=s;i++)
1499    {v=v,L[i][5];}
1500    intvec u=v[1];
1501    for(i=2;i<=s;i++)
1502    {
1503      l=u[size(u)]+v[i];
1504      u=u,l;
1505    }
1506    int m=u[size(u)];
1507    intmat M[m][m];
1508    for(i=1;i<=m;i++)
1509    {
1510      for(j=i+1;j<=m;j++)
1511      {
1512        d=sorting(u,i);
1513        b=sorting(u,j);
1514        if(d==b){k=contact(L[d]);
1515        w1=multsequence(L[d]);
1516        if(size(w1)<k){for(p=size(w1)+1;p<=k;p++)
1517        {w1=w1,1;} }
1518        for(o=1;o<=k;o++){c=c+w1[o]*w1[o];}
1519        M[i,j]=c;
1520        c=0;
1521      }
1522      else
1523      {M[i,j]=intersection(L[d],L[b]);}
1524    }
1525  }
1526  return(M);
1527}
1528/////////////////////////////////////////////////////////////////
1529static proc lengths(matrix M)
1530{
1531  intvec v;
1532  int s,i,j;
1533  for(i=1;i<=nrows(M);i++)
1534  {
1535    for(j=1;j<=ncols(M);j++)
1536    {
1537      if(M[i,j]==t)
1538      {
1539        v[i]=j-1;
1540        if(i==nrows(M)){s=1;}
1541        break;
1542      }
1543    }
1544  }
1545  if(s==0){v[nrows(M)]=-j;}
1546  return(v);
1547}
1548//////////////////////////////////////////////////////////////////////
1549static proc sorting(intvec u,int k)
1550{
1551  int s=size(u);
1552  int i;
1553  for(i=1;i<=s;i++)
1554  {
1555    if(u[i]>=k){break;}
1556  }
1557  return(i);
1558}
1559//////////////////////////////////////////////////////////////////////
1560proc MultiplicitySequence(ideal i)
1561"USAGE:  MultiplicitySequence(i); i ideal
1562ASSUME:  i is the defining ideal of a (reducible) plane curve singularity.
1563RETURN:  list X of charts. Each chart contains the multiplicity sequence of
1564         the corresponding branch.
1565EXAMPLE: example MultiplicitySequence;  shows an example
1566"
1567{
1568  def s=CurveParam(i);
1569  setring s;
1570  int j,k;
1571  def r1=converter(Param);
1572  setring r1;
1573  list Y=hne;
1574  list X;
1575  for(j=1;j<=size(Y);j++)
1576  {
1577    for(k=1;k<=Y[j][5];k++)
1578    {
1579      X=insert(X,multsequence(Y[j]),size(X));
1580    }
1581  }
1582  return(X);
1583}
1584example
1585{
1586  "EXAMPLE:"; echo = 2;
1587  ring r=0,(x,y),ds;
1588  ideal i=x14-x4y7-y11;
1589  MultiplicitySequence(i);
1590}
1591/////////////////////////////////////////////////////////////////////////
1592proc IntersectionMatrix(ideal i)
1593"USAGE:  IntersectionMatrix(i); i ideal
1594ASSUME:  i is the defining ideal of a (reducible) plane curve singularity.
1595RETURN:  intmat of the intersection multiplicities of the branches.
1596EXAMPLE: example IntersectionMatrix;  shows an example
1597"
1598{
1599  def s=CurveParam(i);
1600  setring s;
1601  int j,k;
1602  def r1=converter(Param);
1603  setring r1;
1604  list Y=hne;
1605  return(intermat(Y));
1606}
1607example
1608{
1609  "EXAMPLE:"; echo = 2;
1610  ring r=0,(x,y),ds;
1611  ideal i=x14-x4y7-y11;
1612  IntersectionMatrix(i);
1613}
1614///////////////////////////////////////////////////////////////////////////
1615proc CharacteristicExponents(ideal i)
1616"USAGE:  CharacteristicExponents(i); i ideal
1617ASSUME:  i is the defining ideal of a (reducible) plane curve singularity.
1618RETURN:  list X of charts. Each chart contains the characteristic exponents
1619         of the corresponding branch.
1620EXAMPLE: example CharacteristicExponents;  shows an example
1621"
1622{
1623  def s=CurveParam(i);
1624  setring s;
1625  int j,k;
1626  def r1=converter(Param);
1627  setring r1;
1628  list X;
1629  list Y=hne;
1630  for(j=1;j<=size(Y);j++)
1631  {
1632    for(k=1;k<=Y[j][5];k++)
1633    {
1634      X=insert(X,invariants(Y[j])[1],size(X));
1635    }
1636  }
1637  return(X);
1638}
1639example
1640{
1641  "EXAMPLE:"; echo = 2;
1642  ring r=0,(x,y),ds;
1643  ideal i=x14-x4y7-y11;
1644  CharacteristicExponents(i);
1645}
1646/////////////////////////////////////////////////////////////////////////////
1647static proc contactNumber(int a,intvec v1,intvec v2)
1648{
1649//====  a is the intersection multiplicity of the branches
1650//====  v1,v2 are the multiplicity sequences
1651  int i,c,d;
1652  if(size(v1)>size(v2))
1653  {
1654    for(i=size(v2)+1;i<=size(v1);i++)
1655    {
1656      v2[i]=1;
1657    }
1658  }
1659  if(size(v1)<size(v2))
1660  {
1661    for(i=size(v1)+1;i<=size(v2);i++)
1662    {
1663      v1[i]=1;
1664    }
1665  }
1666  while((c<a)&&(d<size(v1)))
1667  {
1668    d++;
1669    c=c+v1[d]*v2[d];
1670  }
1671  if(c<a)
1672  {
1673    d=d+a-c;
1674  }
1675   return(d);
1676}
1677////////////////////////////////////////////////////////////////////////////
1678proc ContactMatrix(ideal I)
1679"USAGE:  ContactMatrix(I); I ideal
1680ASSUME:  i is the defining ideal of a (reducible) plane curve singularity.
1681RETURN:  intmat N of the contact matrix of the braches of the curve.
1682EXAMPLE: example ContactMatrix;  shows an example
1683"
1684{
1685  def s=CurveParam(I);
1686  setring s;
1687  int j,k,i;
1688  def r1=converter(Param);
1689  setring r1;
1690  list Y=hne;
1691  list L;
1692  for(j=1;j<=size(Y);j++)
1693  {
1694    for(k=1;k<=Y[j][5];k++)
1695    {
1696      L=insert(L,multsequence(Y[j]),size(L));
1697     }
1698  }
1699  intmat M=intermat(Y);
1700  intmat N[nrows(M)][ncols(M)];
1701  for(i=1;i<=nrows(M);i++)
1702  {
1703  for(j=i+1;j<=ncols(M);j++)
1704  {
1705    N[i,j]=contactNumber(M[i,j],L[i],L[j]);
1706    N[j,i]=N[i,j];}
1707  }
1708  return(N);
1709}
1710example
1711{
1712  "EXAMPLE:"; echo = 2;
1713  ring r=0,(x,y),ds;
1714  ideal i=x14-x4y7-y11;
1715  ContactMatrix(i);
1716}
1717///////////////////////////////////////////////////////////////////////////
1718proc plainInvariants(ideal i)
1719"USAGE:  plainInvariants(i); i ideal
1720ASSUME:  i is the defining ideal of a (reducible) plane curve singularity.
1721RETURN:  list L of charts. L[j] is the invariants of the jth branch and the last entry
1722         of L is a list containing the intersection matrix,contact matrix,resolution
1723         graph of the curve.
1724         L[j][1]:     intvec (characteristic exponents of the jth branch)
1725         L[j][2]:    intvec (generators of the semigroup of the jth branch)
1726         L[j][3]:    intvec (first components of the puiseux pairs of the jth branch)
1727         L[j][4]:    intvec (second components of the puiseux pairs of the jth branch)
1728         L[j][5]:    int    (degree of conductor of the jth branch)
1729         L[j][6]:    intvec (multiplicity sequence of the jth branch.
1730         L[last][1]: intmat (intersection matrix of the branches)
1731         L[last][2]: intmat (contact matrix of the branches)
1732         L[last][3]: intmat (resolution graph of the curve)
1733SEE ALSO: MultiplicitySequence, CharacteristicExponents, IntersectionMatrix,
1734          ContactMatrix
1735EXAMPLE: example plainInvariants;  shows an example
1736"
1737{
1738  def s=CurveParam(i);
1739  setring s;
1740  int j,k;
1741  def r1=converter(Param);
1742  setring r1;
1743  list Y=hne;
1744  list L,X,Q;
1745  for(j=1;j<=size(Y);j++)
1746  {
1747    for(k=1;k<=Y[j][5];k++)
1748    {
1749      L=insert(L,invariants(Y[j]),size(L));   //computes the same thing again
1750      X=insert(X,invariants(Y[j])[1],size(X));
1751    }
1752  }
1753  L=insert(L,list(),size(L));
1754  L[size(L)]=insert(L[size(L)],intermat(Y),size(L[size(L)]));
1755  intmat N[nrows(intermat(Y))][ncols(intermat(Y))];
1756  for(k=1;k<=nrows(intermat(Y));k++)
1757  {
1758    for(j=k+1;j<=ncols(intermat(Y));j++)
1759    {
1760      N[k,j]=contactNumber(intermat(Y)[k,j],L[k][6],L[j][6]);
1761      N[j,k]=N[k,j];
1762    }
1763  }
1764  L[size(L)]=insert(L[size(L)],N,size(L[size(L)]));
1765  Q=L[size(L)][2],X;
1766  L[size(L)]=insert(L[size(L)],resolutiongraph(Q),size(L[size(L)]));
1767  return(L);
1768}
1769example
1770{
1771  "EXAMPLE:"; echo = 2;
1772  ring r=0,(x,y),ds;
1773  ideal i=x14-x4y7-y11;
1774  plainInvariants(i);
1775}
1776////////////////////////////////////////////////////////////////////////////////////
1777static proc bound(poly x,poly y)
1778{
1779  poly z=x+y;
1780  int m=ord(z);
1781  int i;
1782  int c=-1;
1783  for(i=2;i<=size(z);i++)
1784  {
1785    if(gcd(m,leadexp(z[i])[1])==1)
1786    {
1787      c=2*leadexp(z[i])[1];
1788      break;
1789    }
1790    else
1791    {
1792      m=gcd(m,leadexp(z[i])[1]);
1793    }
1794  }
1795  return(c);
1796}
1797/////////////////////////////////////////////////////////////////////////////////////
1798static proc boundparam(poly f)
1799{
1800  int i;
1801  int l=size(f);
1802  int m=leadexp(f[l])[1];
1803  for(i=l-1;i>=1;i--)
1804  {
1805    if(gcd(m,leadexp(f[i])[1])==1)
1806    {
1807      i=i-1;
1808      break;
1809    }
1810    else
1811    {
1812      m=gcd(m,leadexp(f[i])[1]);
1813    }
1814  }
1815  return(2*leadexp(f[i+1])[1]);
1816}
Note: See TracBrowser for help on using the repository browser.