source: git/Singular/LIB/hyperel.lib @ 7d56875

fieker-DuValspielwiese
Last change on this file since 7d56875 was 794151, checked in by Hans Schönemann <hannes@…>, 17 years ago
*hannes: curvepar, hyperel git-svn-id: file:///usr/local/Singular/svn/trunk@10078 2c84dea3-7e68-4137-9b89-c4e89433aadc
  • Property mode set to 100644
File size: 24.4 KB
Line 
1///////////////////////////////////////////////////////////////////////////////
2version="$Id: hyperel.lib,v 1.1 2007-05-24 11:57:58 Singular Exp $";
3category="Teaching";
4info="
5LIBRARY:    hyperelliptic.lib
6AUTHOR:     Markus Hochstetter, markushochstetter@gmx.de
7
8 NOTE: This library provides procedures for computing with divisors in the
9       jacobian of hyperelliptic curves. In addition procedures are available
10       for computing the rational representation of divisors and vice versa.
11       The library is intended to be used for teaching and demonstrating
12       purposes but not for efficient computations.
13
14
15
16
17PROCEDURES:
18  ishyper(h,f)           test, if y^2+h(x)y=f(x) is hyperelliptic
19  isoncurve(P,h,f)       test, if point P is on C: y^2+h(x)y=f(x)
20  chinrestp(b,moduli)    compute polynom x, s.t. x=b[i] mod moduli[i]
21  norm(a,b,h,f)          norm of a(x)-b(x)y in IF[C]
22  multi(a,b,c,d,h,f)     (a(x)-b(x)y)*(c(x)-d(x)y)  in IF[C]
23  ratrep (P,h,f)         returns polynomials a,b, s.t. div(a,b)=P
24  divisor(a,b,h,f,[])    computes divisor of a(x)-b(x)y
25  gcddivisor(p,q)        gcd of the divisors p and q
26  semidiv(D,h,f)         semireduced divisor of the pair of polys D[1], D[2]
27  cantoradd(D,Q,h,f)     adding divisors of the hyperell. curve y^2+h(x)y=f(x)
28  cantorred(D,h,f)       returns reduced divisor which is equivalent to D
29  double(D,h,f)          computes 2*D on y^2+h(x)y=f(x)
30  cantormult(m,D,h,f)    computes m*D on y^2+h(x)y=f(x)
31
32              [parameters in square brackets are optional]
33";
34///////////////////////////////////////////////////////////////////////////////
35
36
37//=============== Test, if a given curve is hyperelliptic =====================
38
39proc ishyper(poly h, poly f)
40"USAGE:   ishyper(h,f); h,f poly
41RETURN:  1 if y^2+h(x)y=f(x) is hyperelliptic, else 0
42NOTE:    Tests, if y^2+h(x)y=f(x) is a hyperelliptic curve.
43         Curve is defined over basering. Additionally shows error-messages.
44EXAMPLE: example ishyper; shows an example
45"
46{
47  // constructing a copy of the basering (only variable x),
48  // with variables x,y.
49  def R=basering;
50  list l= ringlist(R);
51  list ll=l[2];
52  ll="x","y";
53  l[2]=ll;
54  intvec v= l[3][1][2];
55  v=v,1;
56  l[3][1][2]=v;
57  def s=ring(l);
58  setring s;
59
60  // test, if y^2 + hy - f is hyperelliptic.
61  int i=1;
62  poly h=imap(R,h);
63  poly f=imap(R,f);
64  poly F=y2 + h*y - f;
65  ideal I=F, diff(F,x) , diff(F,y);
66  ideal J=std(I);
67  if ( J != 1 )
68  {
69      i=0;
70      "The curve is singular!";
71  }
72  if ( deg(f) mod 2 != 1 )
73  {
74     i=0;
75      "The polynomial ",f," has even degree!";
76  }
77  if ( leadcoef(f) != 1 )
78  {
79     i=0;
80     "The polynomial ",f," is not monic!";
81  }
82  if ( 2*deg(h) > deg(f)-1 )
83  {
84    i=0;
85     "The polynomial ",h," has degree ",deg(h),"!";
86  }
87  setring(R);
88  return(i);
89}
90example
91{ "EXAMPLE:"; echo = 2;
92   ring R=7,x,dp;
93   // hyperelliptic curve y^2 + h*y = f
94   poly h=x;
95   poly f=x5+5x4+6x2+x+3;
96   ishyper(h,f);
97}
98
99
100//================= Test, if a given ponit is on the curve ====================
101
102proc isoncurve(list P, poly h, poly f)
103"USAGE:   isoncurve(P,h,f); h,f poly
104RETURN:  1 or 0 (if P is on curve or not)
105NOTE:    Tests, if P=(P[1],P[2]) is on the hyperelliptic curve
106         y^2+h(x)y=f(x) .
107         Curve is defined over basering.
108EXAMPLE: example isoncurve; shows an example
109"
110{
111   if ( P[2]^2 + subst(h,var(1),P[1])*P[2] - subst(f,var(1),P[1]) == 0 )
112   {
113      return(1);
114   }
115   return(0);
116}
117example
118{ "EXAMPLE:"; echo = 2;
119   ring R=7,x,dp;
120   // hyperelliptic curve y^2 + h*y = f
121   poly h=x;
122   poly f=x5+5x4+6x2+x+3;
123   list P=2,3;
124   isoncurve(P,h,f);
125}
126
127
128//====================== Remainder of a polynomial division ===================
129
130proc divrem(f,g)
131"USAGE:   divrem(f,g);  f,g poly
132RETURN:  remainder of the division f/g
133NOTE:    Computes R, s.t. f=a*g + R, and deg(R) < deg(g)
134EXAMPLE: example divrem; shows an example
135"
136{
137   return(reduce(f,std(g)));
138}
139example
140{ "EXAMPLE:"; echo = 2;
141   ring R=0,x,dp;
142   divrem(x2+1,x2);
143}
144
145
146//================ chinese remainder theorem for polynomials ==================
147
148proc chinrestp(list b,list moduli)
149"USAGE:   chinrestp(b,moduli);  b, moduli list of polynomials
150RETURN:  polynomial x, s.t. x= b[i] mod moduli[i]
151NOTE:    chinese remainder theorem for polynomials
152EXAMPLE: example chinrestp; shows an example
153"
154{
155   int i;
156   int n=size(moduli);
157   poly M=1;
158   for(i=1;i<=n;i++)
159   {
160      M=M*moduli[i];
161   }
162   list m;
163   for(i=1;i<=n;i++)
164   {
165      m[i]=M/moduli[i];
166   }
167   list y;
168   for(i=1;i<=n;i++)
169   {
170      y[i]= extgcd(moduli[i],m[i])[3];
171   }
172   poly B=0;
173   for(i=1;i<=n;i++)
174   {
175      B=B+y[i]*m[i]*b[i];
176   }
177   B=divrem(B,M);
178   return(B);
179}
180example
181{ "EXAMPLE:"; echo = 2;
182   ring R=7,x,dp;
183   list b=3x-4, -3x2+1, 1, 4;
184   list moduli=(x-2)^2, (x-5)^3, x-1, x-6;
185   chinrestp(b,moduli);
186}
187
188
189//========================= norm of a polynomial ===============================
190
191proc norm(poly a, poly b, poly h, poly f)
192"USAGE:   norm(a,b,h,f);
193RETURN:  norm of a(x)-b(x)y in IF[C]
194NOTE:    The norm is a polynomial in just one variable.
195         Curve C:y^2+h(x)y=f(x) is defined over basering.
196EXAMPLE: example norm; shows an example
197"
198{
199   poly n=a^2+a*b*h-b^2*f;
200   return(n);
201}
202example
203{ "EXAMPLE:"; echo = 2;
204   ring R=7,x,dp;
205   // hyperelliptic curve y^2 + h*y = f
206   poly h=x;
207   poly f=x5+5x4+6x2+x+3;
208   poly a=x2+1;
209   poly b=x;
210   norm(a,b,h,f);
211}
212
213
214//========== multiplikation of polynomials in the coordinate ring =============
215
216proc multi(poly a, poly b, poly c, poly d, poly h, poly f)
217"USAGE:   multi(a,b,c,d,h,f);
218RETURN:  list L with L[1]-L[2]y=(a(x)-b(x)y)*(c(x)-d(x)y) in IF[C]
219NOTE:    Curve C:y^2+h(x)y=f(x) is defined over basering.
220EXAMPLE: example multi; shows an example
221"
222{
223   poly A=a*c + b*d*f;
224   poly B=b*c +a*d + b*h*d;
225   return (list(A,B));
226}
227example
228{ "EXAMPLE:"; echo = 2;
229   ring R=7,x,dp;
230   poly h=x;
231   poly f=x5+5x4+6x2+x+3;
232   // hyperelliptic curve y^2 + h*y = f
233   poly a=x2+1;
234   poly b=x;
235   poly c=5;
236   poly d=-x;
237   multi(a,b,c,d,h,f);
238}
239
240
241//================== polynomial expansion around a point ========================
242
243proc darst(list P,int k, poly h, poly f)
244"USAGE:   darst(P,k,h,f);
245RETURN:  list c of length k
246NOTE:    expansion around point P in IF[C], s.t.
247         y=c[1]+c[2]*(x-P[1]) +...+c[k]*(x-P[1])^k-1 + rest.
248         Curve C:y^2+h(x)y=f(x) is defined over basering.
249EXAMPLE: example darst; shows an example
250"
251{
252
253   if ( P[2] == -P[2]- subst(h,var(1),P[1]))
254   {
255    ERROR("no special points allowed");
256   }
257   list c;
258   list r;
259   list n;
260   poly N;
261   c[1]=P[2];
262   r[1]=list(0,-1,1,0);
263   poly r1,r2,r3,r4;
264   // rational function are represented as (r1 - r2*y) / (r3 - r4*y)
265
266   for (int i=1; i<k ; i++)
267   {
268      r1=r[i][1]-c[i]*r[i][3];
269      r2=r[i][2]-c[i]*r[i][4];
270      r3=r[i][3];
271      r4=r[i][4];
272      n=multi(r3,r4,r1+r2*h,-r2,h,f);
273      N=r1*r1 + r1*r2*h-r2*r2*f;
274      r[i+1]=list(N/(var(1)-P[1]),0,n[1],n[2]);
275      while ((divrem(r[i+1][1],var(1)-P[1]) ==0) and (divrem(r[i+1][2],var(1)-P[1]) ==0) and (divrem(r[i+1][3],var(1)-P[1]) ==0) and (divrem(r[i+1][4],var(1)-P[1]) ==0))
276      {
277          // reducing the rationl function
278          //(r[i+1][1] - r[i+1][2]*y)/(r[i+1][3] - r[i+1][4]) , otherwise there
279          // could be a pole, s.t. conditions are not fulfilled.
280          r[i+1][1]=(r[i+1][1]) / (var(1)-P[1]);
281          r[i+1][2]=(r[i+1][2]) / (var(1)-P[1]);
282          r[i+1][3]=(r[i+1][3]) / (var(1)-P[1]);
283          r[i+1][4]=(r[i+1][4]) / (var(1)-P[1]);
284      }
285      c[i+1]=(subst(r[i+1][1],var(1),P[1]) - subst(r[i+1][2],var(1),P[1])*P[2]) / (subst(r[i+1][3],var(1),P[1]) - subst(r[i+1][4],var(1),P[1])*P[2]);
286   }
287   return(c);
288}
289example
290{ "EXAMPLE:"; echo = 2;
291   ring R=7,x,dp;
292   // hyperelliptic curve y^2 + h*y = f
293   poly h=x;
294   poly f=x5+5x4+6x2+x+3;
295   list P=5,3;
296   darst(P,3,h,f);
297}
298
299
300//================ rational representation of a divisor =======================
301
302proc ratrep1 (list P, poly h, poly f)
303"USAGE:   ratrep1(P,k,h,f);
304RETURN:  list (a,b)
305NOTE:    Important: P has to be semireduced!
306         Computes rational representation of the divisor
307         P[1][3]*(P[1][1], P[1][2]) +...+ P[sizeof(P)][3]*
308         *(P[sizeof(P)][1], P[sizeof(P)][2]) - (*)infty=div(a,b)
309         Divisor P has to be semireduced.
310         Curve C:y^2+h(x)y=f(x) is defined over basering.
311SEE AlSO: ratrep
312EXAMPLE: example ratrep1; shows an example
313"
314{
315   poly a=1;
316   list b;
317   list m;
318   list koef;
319   int k;
320
321   // Determination of the polynomial b[i] for each point using procedure darst
322   for (int i=1 ; i<= size(P); i++)
323   {
324      a=a*(var(1)-P[i][1])^(P[i][3]);    // computing polynomial a
325      m[i]=(var(1)-P[i][1])^(P[i][3]);
326      b[i]=P[i][2];
327      k=1;
328
329      while (divrem(b[i]*b[i] + b[i] *h - f,(var(1)-P[i][1])^(P[i][3])) != 0)
330      {
331          k=k+1;   // b[i]=P[i][2];
332          koef=darst(list (P[i][1],P[i][2]), k, h,f);
333          // could be improved, if one doesn't compute list coef completely new
334          // every time
335          b[i]=b[i]+ koef[k]*(var(1)-P[i][1])^(k-1);
336      }
337   }
338   // Return polynomial a and b. Polynomial b is solution of the congruencies
339   // b[i] mod m[i] .
340   return(list(a,chinrestp(b,m)));
341
342}
343example
344{ "EXAMPLE:"; echo = 2;
345   ring R=7,x,dp;
346   // hyperelliptic curve y^2 + h*y = f
347   poly h=x;
348   poly f=x5+5x4+6x2+x+3;
349   //divisor P
350   list P=list(-1,-3,1),list(1,1,1);
351   ratrep1(P,h,f);
352}
353
354
355//================ rational representation of a divisor =======================
356
357proc ratrep (list P, poly h, poly f)
358"USAGE:   ratrep(P,k,h,f);
359RETURN:  list (a,b)
360NOTE:    Importatnt: P has to be semireduced!
361         Computes rational representation of the divisor
362         P[1][3]*(P[1][1], P[1][2]) +...+ P[sizeof(P)][3]*
363         *(P[sizeof(P)][1], P[sizeof(P)][2]) - (*)infty=div(a,b)
364         Divisor P has to be semireduced.
365         Curve C:y^2+h(x)y=f(x) is defined over basering.
366         Works faster than ratrep1.
367SEE ALSO: ratrep1
368EXAMPLE: example ratrep; shows an example
369"
370{
371   poly a=1;
372   list b;
373   list m;
374   list koef;
375   int k;
376
377   poly c;
378   list r;
379   list n;
380   poly Norm;
381   poly r1,r2,r3,r4;
382
383   // Determination of the polynomial b[i] for each point using procedure darst
384   for (int i=1 ; i<= size(P); i++)
385   {
386      a=a*(var(1)-P[i][1])^(P[i][3]);     // computing polynomial a
387      m[i]=(var(1)-P[i][1])^(P[i][3]);
388      b[i]=P[i][2];
389      k=1;
390      c=P[i][2];
391      r=0,-1,1,0;
392      while (divrem(b[i]*b[i] + b[i] *h - f,(var(1)-P[i][1])^(P[i][3])) != 0)
393      {
394          k=k+1;
395          // here, the procedure darst was integrateg. In every pass a new
396          // coefficient c[i] is determined.
397          r1=r[1]-c*r[3];
398          r2=r[2]-c*r[4];
399          r3=r[3];
400          r4=r[4];
401          n=multi(r3,r4,r1+r2*h,-r2,h,f);
402          Norm=r1*r1 + r1*r2*h-r2*r2*f;
403          r=list(Norm/(var(1)-P[i][1]),0,n[1],n[2]);
404          while ((divrem(r[1],var(1)-P[i][1]) ==0) and (divrem(r[2],var(1)-P[i][1]) ==0) and (divrem(r[3],var(1)-P[i][1]) ==0) and (divrem(r[4],var(1)-P[i][1]) ==0))
405          {
406
407           // reducing the rationl function
408           // (r[1]-r[2]y)/(r[3]-r[4]y) , otherwise there
409           // could be a pole, s.t. conditions are not fulfilled.
410              r[1]=(r[1]) / (var(1)-P[i][1]);
411              r[2]=(r[2]) / (var(1)-P[i][1]);
412              r[3]=(r[3]) / (var(1)-P[i][1]);
413              r[4]=(r[4]) / (var(1)-P[i][1]);
414          }
415          c=(subst(r[1],var(1),P[i][1]) - subst(r[2],var(1),P[i][1])*P[i][2]) / (subst(r[3],var(1),P[i][1]) - subst(r[4],var(1),P[i][1])*P[i][2]);
416          b[i]=b[i]+ c*(var(1)-P[i][1])^(k-1);
417      }
418   }
419   // return polynomial a and b. Polynomial b is solution of the congruencies
420   // b[i] mod m[i] .
421   return(list(a,chinrestp(b,m)));
422}
423example
424{ "EXAMPLE:"; echo = 2;
425   ring R=7,x,dp;
426   // hyperelliptic curve y^2 + h*y = f
427   poly h=x;
428   poly f=x5+5x4+6x2+x+3;
429   //Divisor P
430   list P=list(-1,-3,1),list(1,1,1);
431   ratrep(P,h,f);
432}
433
434
435//============== Order of a zero in a polynomial ==============================
436
437proc ordnung(poly x0 , poly g)
438"USAGE:   ordnung(x0,g);
439RETURN:  int i
440NOTE:    i is maximal, s.t. (x-x0)^i divides g
441EXAMPLE: example ordnung; shows an example
442"
443{
444   poly gg=g;
445   int i;
446   while ( divrem(gg,var(1)-x0) ==0 )
447   {
448      i=i+1;
449      gg=gg/(var(1)-x0);
450   }
451   return(i);
452}
453example
454{ "EXAMPLE:"; echo = 2;
455   ring R=0,x,dp;
456   poly g=(x-5)^7*(x-3)^2;
457   number x0=5;
458   ordnung(x0,g);
459}
460
461
462//================== divisor of a polynomial function =========================
463
464proc divisor(poly a, poly b, poly h, poly f, list #)
465"USAGE:   divisor(a,b,h,f); optional: divisor(a,b,h,f,s); s=0,1
466RETURN:  list P
467NOTE:    P[1][3]*(P[1][1], P[1][2]) +...+ P[sizeof(P)][3]*
468         *(P[sizeof(P)][1], P[sizeof(P)][2]) - (*)infty=div(a(x)-b(x)y)
469         if there is an optional parameter s!=0, then divisor additonally
470         returns a parameter, which says, whether irreducible polynomials
471         occured during computations or not. Otherwise only warnings are
472         displayed on the monitor. For s=0 nothing happens.
473         Curve C:y^2+h(x)y=f(x) is defined over basering.
474EXAMPLE: example divisor; shows an example
475"
476{
477   list p;
478   int j;
479   poly x0;
480   list y;
481   list fa=factorize(gcd(a,b));  // wanted: common roots of a and b
482   poly Norm=norm(a,b,h,f);
483
484   int s;
485   int irred=0;
486   if (size(#)>0)
487   {
488      s=#[1];
489   }
490   else
491   {
492      s=0;
493   }
494
495   for (int i=2; i<=size(fa[1]) ; i++)
496   {
497      // searching roots by finding polynomials of degree 1
498      if ( deg(fa[1][i]) !=1 )
499      {
500         if (s==0)
501         {
502         "WARNIG: ", fa[1][i], "is irreducible over this field !";
503         }
504         else
505         {
506            irred=1;
507         }
508      }
509      else
510      {
511         x0=var(1) - fa[1][i];
512         // finding the y-coordinates; max. 2
513         y= factorize(var(1)^2 + var(1)*subst(h,var(1),x0) - subst(f,var(1),x0));
514         if ( deg(y[1][2]) == 1)
515         // if root belongs to point on curve, then...
516         {
517            // compute order of a-b*y in the founded point
518            j=j+1;
519            p[j]=list(x0,var(1)-y[1][2],fa[2][i]);
520            if ( y[2][2]== 1)          // ordinary point
521            {
522               j=j+1;
523               p[j]=list(x0 , var(1)-y[1][3] , fa[2][i] );
524               if (a/(var(1)-x0)^(fa[2][i]) - b/(var(1)-x0)^(fa[2][i]) * p[j][2] ==0 )
525               {
526                  p[j][3]= p[j][3] + ordnung(x0,norm(a/(var(1)-x0)^(fa[2][i]) , b/(var(1)-x0)^(fa[2][i]),h,f));
527               }
528               if (a/(var(1)-x0)^(fa[2][i]) - b/(var(1)-x0)^(fa[2][i]) * p[j-1][2] ==0 )
529               {
530               p[j-1][3]=p[j-1][3] + ordnung(x0,norm(a/(var(1)-x0)^(fa[2][i]) , b/(var(1)-x0)^(fa[2][i]),h,f));
531               }
532            }
533            else        // special point
534            {
535               p[j][3]=p[j][3] *2  ;
536               if (a/(var(1)-x0)^(fa[2][i]) - b/(var(1)-x0)^(fa[2][i]) * p[j][2] ==0 )
537               {
538                  p[j][3]=p[j][3] + ordnung(x0,norm(a/(var(1)-x0)^(fa[2][i]) , b/(var(1)-x0)^(fa[2][i]),h,f));
539               }
540            }
541
542         }
543         // Norm of a-b*y is reduced by common root of a and b
544         // (is worked off)
545         Norm = Norm/((var(1)-x0)^(ordnung(x0,Norm)));
546      }
547   }
548
549   // some points are still missing; points for which a and b have no common
550   // roots, but norm(a-b*Y)=0 .
551   fa=factorize(Norm);
552   for ( i=2 ; i<=size(fa[1]) ; i++)
553   {
554      if ( deg(fa[1][i]) !=1)
555      {
556         if (s==0)
557         {
558         "WARNING: ", fa[1][i], "is irreducible over this field !";
559         }
560         else
561         {
562            irred=1;
563         }
564      }
565      else
566      {
567         x0=var(1) - fa[1][i];
568         y= factorize(var(1)^2 + var(1)*subst(h,var(1),x0) - subst(f,var(1),x0));
569         if ( deg(y[1][2]) == 1)
570         // if root belongs to point on curve, then...
571         {
572            if (subst(a,var(1),x0)- subst(b,var(1),x0)* (var(1)-y[1][2]) ==0)
573            {
574               p[size(p)+1]=list(x0,var(1)-y[1][2], ordnung(x0,Norm,h,f));
575            }
576            if ( y[2][2]== 1)     // ordinary point
577            {
578               if (subst(a,var(1),x0)- subst(b,var(1),x0)* (var(1)-y[1][3]) ==0)
579               {
580               p[size(p)+1]=list(x0 , var(1)-y[1][3] , ordnung(x0,Norm,h,f));
581               }
582            }
583         }
584      }
585   }
586   if (s==0)
587   {
588      return(p);
589   }
590   return(p,irred);
591}
592example
593{ "EXAMPLE:"; echo = 2;
594   ring R=7,x,dp;
595   // hyperelliptic curve y^2 + h*y = f
596   poly h=x;
597   poly f=x5+5x4+6x2+x+3;
598   poly a=(x-1)^2*(x-6);
599   poly b=0;
600   divisor(a,b,h,f,1);
601}
602
603
604//===================== gcd of two divisors ===================================
605
606proc gcddivisor(list p, list q)
607"USAGE:   gcddivisor(p,q);
608RETURN:  list P
609NOTE:    gcd of two divisors
610EXAMPLE: example gcddivisor; shows an example
611"
612{
613   list e;
614   int i,j;
615   for (i=1 ; i<= size(p) ; i++)
616   {
617      for (j=1 ; j<= size(q) ; j++)
618      {
619         if ( p[i][1] == q[j][1] and p[i][2] == q[j][2])
620         {
621            if ( p[i][3] <= q[j][3] )
622            {
623               e[size(e)+1]= list (p[i][1] , p[i][2] , p[i][3]);
624            }
625            else
626            {
627               e[size(e)+1]= list (q[j][1] , q[j][2] , q[j][3]);
628            }
629         }
630      }
631   }
632   return(e);
633}
634example
635{ "EXAMPLE:"; echo = 2;
636   ring R=7,x,dp;
637   // hyperelliptic curve y^2 + h*y = f
638   poly h=x;
639   poly f=x5+5x4+6x2+x+3;
640   // two divisors
641   list p=list(-1,-3,1),list(1,1,2);
642   list q=list(1,1,1),list(2,2,1);
643   gcddivisor(p,q);
644}
645
646
647//========== semireduced divisor from rational representation=================
648
649proc semidiv(list D,poly h, poly f)
650"USAGE:   semidiv(D,h,f);
651RETURN:  list P
652NOTE:    important: Divisor D has to be semireduced!
653         Computes semireduced divisor
654         P[1][3]*(P[1][1], P[1][2]) +...+ P[sizeof(P)][3]*
655         *(P[sizeof(P)][1], P[sizeof(P)][2]) - (*)infty=div(D[1],D[2])
656         Curve C:y^2+h(x)y=f(x) is defined over basering.
657EXAMPLE: example semidiv; shows an example
658"
659{
660   if ( deg(D[2]) >= deg(D[1]) or divrem(D[2]^2+D[2]*h-f,D[1]) != 0 )
661   {
662      ERROR("Pair of polynomials doesn't belong to semireduced divisor!");
663   }
664   list D1,D2;
665   int s1,s2;
666   D1,s1=divisor(D[1],0,h,f,1);
667   D2,s2=divisor(D[2],1,h,f,1);
668
669   // Only if irreducible polynomials occured in D1 !and! D2 a warning
670   // is necessary.
671   if (s1==1 and s2==1)
672   {
673      "Attention:
674     Perhaps some points were not found over this field!";
675   }
676   return(gcddivisor(D1,D2));
677}
678example
679{ "EXAMPLE:"; echo = 2;
680   ring R=7,x,dp;
681   // hyperelliptic curve y^2 + h*y = f
682   poly h=x;
683   poly f=x5+5x4+6x2+x+3;
684   // Divisor
685   list D=x2-1,2x-1;
686   semidiv(D,h,f);
687}
688
689
690//=============== Cantor's algorithm - composition ============================
691
692proc cantoradd(list D, list Q, poly h, poly f)
693"USAGE:   cantoradd(D,Q,h,f);
694RETURN:  list P
695NOTE:    Cantor's Algorithm - composition
696         important: D and Q have to be semireduced!
697         Computes semireduced divisor
698         div(P[1],P[2])= div(D[1],D[2]) + div(Q[1],Q[2])
699         The divisors are defined over the basering.
700         Curve C:y^2+h(x)y=f(x) is defined over the basering.
701EXAMPLE: example cantoradd; shows an example
702"
703{
704   poly a;
705   poly b;
706   list e=extgcd(D[1],Q[1]);
707   if ( e[1]==1 )
708   {
709      a=D[1]*Q[1];
710      b=divrem( e[2]*D[1]*Q[2] + e[3]*Q[1]*D[2] ,a);
711      return(list(a,b));
712   }
713   list c=extgcd(e[1],D[2]+Q[2]+h);
714   poly s1=e[2]*c[2];
715   poly s2=c[2]*e[3];
716   poly s3=c[3];
717   a=D[1]*Q[1]/c[1]^2;
718   b=divrem((s1*D[1]*Q[2] + s2*Q[1]*D[2] + s3*(D[2]*Q[2] + f))/c[1],a);
719   return(list(a,b));
720}
721example
722{ "EXAMPLE:"; echo = 2;
723   ring R=7,x,dp;
724   // hyperelliptic curve y^2 + h*y = f
725   poly h=x;
726   poly f=x5+5x4+6x2+x+3;
727   // two divisors in rational representation
728   list D=x2-1,2x-1;
729   list Q=x2-3x+2,-3x+1;
730   cantoradd(D,Q,h,f);
731}
732
733
734//==================== Cantor's algorithm - reduction =========================
735
736proc cantorred(list D,poly h,poly f)
737"USAGE:   cantorred(D,h,f);
738RETURN:  list N
739NOTE:    Cantor's algorithm - reduction.
740         important: Divisor D has to be semireduced!
741         Computes reduced divisor div(N[1],N[2])= div(D[1],D[2]).
742         The divisors are defined over the basering.
743         Curve C:y^2+h(x)y=f(x) is defined over the basering.
744EXAMPLE: example cantorred; shows an example
745"
746{
747   list N=D;
748   while ( 2*deg(N[1]) > deg(f)-1 )
749   {
750      N[1]=(f - N[2]*h - N[2]^2)/N[1];
751      N[2]=divrem(-h-N[2],N[1]);
752   }
753   N[1]=N[1]/leadcoef(N[1]);
754   return(N);
755}
756example
757{ "EXAMPLE:"; echo = 2;
758   ring R=7,x,dp;
759   // hyperelliptic curve y^2 + h*y = f
760   poly h=x;
761   poly f=x5+5x4+6x2+x+3;
762   // semireduced divisor
763   list D=2x4+3x3-3x-2, -x3-2x2+3x+1;
764   cantorred(D,h,f);
765}
766
767
768//================= doubling a semireduced divisor ============================
769
770proc double(list D, poly h, poly f)
771"USAGE:   double(D,h,f);
772RETURN:  list Q=2*D
773NOTE:    important: Divisor D has to be semireduced!
774         Special case of Cantor's algorithm.
775         Computes reduced divisor div(Q[1],Q[2])= 2*div(D[1],D[2])
776         The divisors are defined over the basering.
777         Curve C:y^2+h(x)y=f(x) is defined over the basering.
778EXAMPLE: example double; shows an example
779"
780{
781   list c=extgcd(D[1], 2*D[2] + h);
782   poly a=D[1]*D[1]/c[1]^2;
783   poly b=divrem((c[2]*D[1]*D[2] + c[3]*(D[2]*D[2] + f))/c[1],a);
784   return(cantorred(list(a,b),h,f));
785}
786example
787{ "EXAMPLE:"; echo = 2;
788   ring R=7,x,dp;
789   // hyperelliptic curve y^2 + h*y = f
790   poly h=x;
791   poly f=x5+5x4+6x2+x+3;
792   // reduced divisor
793   list D=x2-1,2x-1;
794   double(D,h,f);
795}
796
797
798//================ multiples of a semireduced divisor =========================
799
800proc cantormult(int m, list D, poly h, poly f)
801"USAGE:   cantormult(m,D,h,f);
802RETURN:  list res=m*D
803NOTE:    important: Divisor D has to be semireduced!
804         Uses repeated doublings for a faster computation
805         of the reduced divisor m*D.
806         Attention: Factor m is int, this means bounded.
807         For m<0 the inverse of m*D is returned.
808         The divisors are defined over the basering.
809         Curve C:y^2+h(x)y=f(x) is defined over the basering.
810EXAMPLE: example cantormult; shows an example
811"
812{
813   list res=1,0;
814   list bas=D;
815   int exp=m;
816   if (exp==0) { return(list(1,0)); }
817   if (exp==1) { return(D); }
818   if (exp==-1) { return(list(D[1],-D[2]-h)) ; }
819   if ( exp < 0)
820   {
821      exp=-exp;
822   }
823   while ( exp > 0 )
824   {
825     if ( (exp mod 2) !=0 )
826     {
827       res = cantorred(cantoradd(res,bas,h,f),h,f);
828       exp=exp-1;
829     }
830     bas=double(bas,h,f);
831     exp=exp/2;
832   }
833   if ( m < 0 )
834   {
835      res[2]=-res[2]-h;
836   }
837   return(res);
838}
839example
840{ "EXAMPLE:"; echo = 2;
841   ring R=7,x,dp;
842   // hyperelliptic curve y^2 + h*y = f
843   poly h=x;
844   poly f=x5+5x4+6x2+x+3;
845   // reduced divisor
846   list D=x2-1,2x-1;
847   cantormult(34,D,h,f);
848}
849
850
851/*
852//=============================================================================
853//   In the following you find a large example, which demostrates the use of
854//   the most important procedures.
855//=============================================================================
856//---- field with 2^5=32 elements ----
857ring r=(2,a),x,dp;
858minpoly=a5+a2+1;
859
860//---- hyperelliptic curve y^2 + hy = f ----
861poly h=x2+x;
862poly f=x5+x3+1;
863
864//---- two divisors ----
865list l1=list(a30,0,1),list(0,1,1);
866list l2=list(a30,a16,1),list(1,1,1);
867
868//---- their rational representation ----
869list D1=ratrep(l1,h,f); D1;
870//[1]:
871//   x2+(a4+a)*x
872//[2]:
873//   (a)*x+1
874
875list D2=ratrep(l2,h,f); D2;
876//[1]:
877//   x2+(a4+a+1)*x+(a4+a)
878//[2]:
879//   (a3+a2+a+1)*x+(a3+a2+a)
880
881//---- back to the point-based-representation ----
882semidiv(D1,h,f);
883//[1]:
884//   [1]:
885//      (a4+a)
886//   [2]:
887//      0
888//   [3]:
889//      1
890//[2]:
891//   [1]:
892//      0
893//   [2]:
894//      1
895//   [3]:
896//      1
897
898semidiv(D2,h,f);
899//[1]:
900//   [1]:
901//      (a4+a)
902//   [2]:
903//      (a4+a3+a+1)
904//   [3]:
905//      1
906//[2]:
907//   [1]:
908//      1
909//   [2]:
910//      1
911//   [3]:
912//      1
913
914//---- adding D1 and D2 ----
915list D12=cantorred(cantoradd(D1,D2,h,f),h,f);  D12;
916//[1]:
917//   x2+x
918//[2]:
919//   1
920
921//---- D1+D2 in point-based-representation ----
922semidiv(D12,h,f);
923//[1]:
924//   [1]:
925//      1
926//   [2]:
927//      1
928//   [3]:
929//      1
930//[2]:
931//   [1]:
932//      0
933//   [2]:
934//      1
935//   [3]:
936//      1
937
938//---- D1 + D1   (2 possible ways, same result) ----
939cantorred(cantoradd(D1,D1,h,f),h,f);
940double(D1,h,f);
941//[1]:
942//   x2+(a3+1)
943//[2]:
944//   (a4+a3+a+1)*x+(a4+a3+a2+a+1)
945
946//---- order of D1 in the jacobian over the basering ----
947int i=1;
948list E=D1;
949while (E[1] != 1 or E[2] != 0 )
950{
951   E= cantorred(cantoradd(E,D1,h,f),h,f);
952   i=i+1;
953}
954i;   // 482
955
956//---- proof with multiplikation validates the result ----
957cantormult(i,D1,h,f);
958//[1]:
959//   1
960//[2]:
961//   0
962
963//---- computing the inverse of D1 ----
964list d1= cantormult(-1,D1,h,f);  d1;
965//[1]:
966//   x2+(a4+a)*x
967//[2]:
968//   x2+(a+1)*x+1
969
970//---- proof validates the result ----
971cantoradd(d1,D1,h,f);
972//[1]:
973//   1
974//[2]:
975//   0
976
977
978*/
Note: See TracBrowser for help on using the repository browser.