source: git/Singular/LIB/paraplanecurves.lib @ 1628b14

spielwiese
Last change on this file since 1628b14 was 1628b14, checked in by Frank Seelisch <seelisch@…>, 13 years ago
included outputs for time consumption git-svn-id: file:///usr/local/Singular/svn/trunk@14044 2c84dea3-7e68-4137-9b89-c4e89433aadc
  • Property mode set to 100644
File size: 100.5 KB
Line 
1///////////////////////////////////////////////////////////////////////////////
2version="$Id$";
3category="Algebraic Geometry";
4info="
5LIBRARY:  paraplanecurves.lib Rational parametrization of rational plane curves
6
7AUTHORS:  J. Boehm, j.boehm at mx.uni-saarland.de @*
8          W. Decker, decker at mathematik.uni-kl.de> @*
9          S. Laplagne, slaplagn at dm.uba.ar @*
10          F. Seelisch, seelisch at mathematik.uni-kl.de
11
12OVERVIEW:
13
14Suppose C = {f(x,y,z)=0} is a rational plane curve, where f is homogeneous
15of degree n with coefficients in Q and absolutely irreducible (these
16conditions are checked automatically.) @*
17After a first step, realized by a projective automorphism in the procedure
18adjointIdeal, C satisfies: @*
19- C does not have singularities at infinity z=0. @*
20- C does not contain the point (0:1:0) (that is, the dehomogenization of f
21  with respect to z is monic as a polynomial in y). @*
22Considering C in the chart z<>0, the algorithm regards x as transcendental
23and y as algebraic and computes an integral basis in C(x)[y] of the integral
24closure of C[x] in C(x,y) using the  normalization algorithm from
25@ref{normal_lib}: see @ref{integralbasis_lib}. In a future edition of the
26library, also van Hoeij's algorithm for computing the integral basis will
27be available. @*
28From the integral basis, the adjoint ideal is obtained by linear algebra.
29Alternatively, the algorithm starts with a local analysis of the singular
30locus of C. Then, for each  primary component of the singular locus which
31does not correspond to ordinary multiple points or cusps, the integral
32basis algorithm is applied separately. The ordinary multiple points and
33cusps, in turn, are addressed by a straightforward direct algorithm. The
34adjoint ideal is obtained by intersecting all ideals obtained locally.
35The local variant of the algorithm is used by default. @*
36The linear system corresponding to the adjoint ideal maps the curve
37birationally to a rational normal curve in P^(n-2). @*
38Iterating the anticanonical map, the algorithm projects the rational normal
39curve to PP1 for n odd resp. to a conic C2 in PP2 for n even. @*
40In case n is even, the algorithm tests whether there is a rational point on
41C2 and if so gives a parametrization of C2 which is defined over Q. Otherwise,
42the parametrization given is defined over a quadratic field extension of Q. @*
43By inverting the birational map of C to PP1 resp. to C2, a parametrization
44of C is obtained (defined over Q or the quadratic field extension).
45
46REFERENCES:
47
48Janko Boehm: Parametrisierung rationaler Kurven, Diploma Thesis,
49http://www.math.uni-sb.de/ag/schreyer/jb/diplom%20janko%20boehm.pdf
50
51Theo de Jong: An algorithm for computing the integral closure,
52Journal of Symbolic Computation 26 (3) (1998), p. 273-277
53
54Gert-Martin Greuel, Santiago Laplagne, Frank Seelisch: Normalization of Rings,
55Journal of Symbolic Computation 9 (2010), p. 887-901
56
57Mark van Hoeij: An Algorithm for Computing an Integral Basis in an Algebraic
58Function Field, Journal of Symbolic Computation 18 (1994), p. 353-363,
59http://www.math.fsu.edu/~hoeij/papers/comments/jsc1994.html
60
61KEYWORDS:
62Curves; Parametrization; Rational curves; Adjoint ideal; Geometric genus
63
64
65PROCEDURES:
66
67adjointIdeal(poly, [...]);        Adjoint ideal of a plane curve
68invertBirMap(ideal,ideal);        Invert a birational map of algebraic
69                                  varieties
70paraPlaneCurve(poly, [...]);      Compute a rational parametrization of a
71                                  rational plane curve
72rncAntiCanonicalMap(ideal);       Anticanonical map of a rational normal curve
73rationalPointConic(poly);         Finds a point on the conic. This point has
74                                  either coefficients in Q or in a quadratic
75                                  extension field of Q
76mapToRatNormCurve(poly,ideal);    Map a plane rational curve to a rational
77                                  normal curve (RNC)
78rncItProjOdd(ideal);              Map a RNC via successive anticanonical maps
79                                  to PP1
80rncItProjEven(ideal);             Map a RNC via successive anticanonical maps
81                                  to a conic in PP2
82paraConic(poly);                  Compute a rational parametrization of a conic
83testParametrization(poly,ring);   Checks whether a given curve is parametrized
84                                  by a given rational map (defined in the
85                                  given ring)
86testPointConic(poly,ring);        Checks whether a given point (defined in the
87                                  given ring) lies on the given conic.
88";
89
90LIB "elim.lib";
91LIB "general.lib";
92LIB "primdec.lib";
93LIB "absfact.lib";
94LIB "matrix.lib";
95LIB "random.lib";
96LIB "homolog.lib";
97LIB "integralbasis.lib";
98LIB "normal.lib";
99
100
101///////////////////////////////////////////////////////////////////////////////
102proc invertBirMap(ideal phi, ideal I)
103"USAGE: invertBirMap(phi, I); phi ideal, I ideal
104ASSUME: The ideal phi in the basering R represents a birational map of the
105        variety given by the ideal I in R to its image in projective space
106        P = PP^(size(phi)-1).
107NOTE:   The procedure might fail or give a wrong output if phi does
108        not define a birational map.
109RETURN: ring, the coordinate ring of P, with an ideal named J and an ideal
110        named psi.@*
111        The ideal J defines the image of phi.@*
112        The ideal psi gives the inverse of phi.@*
113        Note that the entries of psi should be considered as representatives
114        of classes in the quotient ring R/J.@*
115THEORY: We compute the ideal I(G) in R**S of the graph G of phi.@*
116        The ideal J is given by the intersection of I(G) with S.@*
117        The map psi is given by a relation mod J of those relations
118        in I(G) which are linear in the variables of R.@*
119KEYWORDS: birational map, image, inverse.
120EXAMPLE: example invertBirMap; shows an example
121"
122{
123   def Roriginal = basering;
124   int n = nvars(Roriginal);
125   int m = size(phi);
126   /*phi: P^(n-1) --> P^(m-1)*/
127   list rl = ringlist(Roriginal);
128   int k;
129   for(k = 1; k <= n; k++)
130     {
131        rl[2][k] = "x("+string(k)+")";
132     }
133   for(k = 1; k <= m; k++)
134     {
135        rl[2][k+n] = "y("+string(k)+")";
136     }
137   rl[3]= list(list("dp",1:(n+m)),list("C",0));
138   /*Use Hilbert driven Buchberger*/
139   def Rbig0 = ring(rl);
140   setring Rbig0;
141   ideal I = fetch(Roriginal,I);
142   ideal phi = fetch(Roriginal,phi);
143   ideal mi = maxideal(1);
144   ideal xv = mi[1..n];
145   ideal yv  = mi[n+1..n+m];
146   matrix HM[2][m] = concat(transpose(yv),transpose(phi));
147   ideal graph = sat(I+minor(HM,2),phi)[1];
148   graph = sat(graph,xv)[1];
149   intvec Hgraph = hilb(graph,1);
150   setring Roriginal;
151   rl[3]= list(list("dp",1:n),list("dp",1:m),list("C",0));
152   def Rbig = ring(rl);
153   setring Rbig;
154   ideal graph = imap(Rbig0,graph);
155   graph = std(graph,Hgraph);
156   ideal xv = imap(Rbig0,xv);
157   /*The ideal J defines the image of phi*/
158   ideal J = graph;
159   for(k = 1; k <= n; k++)
160     {
161        J = subst(J,xv[k],0);
162     }
163   J = compress(J);
164   /*now we start inverting phi to psi*/
165   matrix relpsi = diff(xv,graph);
166   for(k = 1; k <= n; k++)
167     {
168        relpsi = subst(relpsi,xv[k],0);
169     }
170    relpsi = compress(relpsi);
171    list rl = ringlist(Rbig);
172    list rl2 = rl[2];
173    rl[2] = list(rl2[n+1..n+m]);
174    rl[3]= list(list("dp",1:m),list("C",0));
175    def Rtarget = ring(rl);
176    setring Rtarget;
177    ideal J = imap(Rbig,J);
178    qring RtargetmodJ = std(J);
179    matrix relpsi = imap(Rbig,relpsi);
180    relpsi = syz(transpose(relpsi));
181    ideal psi = submat(relpsi,1..nrows(relpsi),1);
182    setring Rtarget;
183    ideal psi = imap(RtargetmodJ,psi);
184    export(J,psi);
185    int p = printlevel - voice + 3;
186    dbprint(p,"// 'invertBirMap' created a ring together with two ideals J and psi.");
187    dbprint(p,"// Supposing you typed, say,  def RPn = invertBirMap(phi,I);");
188    dbprint(p,"// you may access the ideals by typing");
189    dbprint(p,"//      setring RPn; J; psi;");
190    return(Rtarget);
191 }
192
193example
194{ "EXAMPLE:"; echo=2;
195  ring R = 0,(x,y,z),dp;
196  poly f = y^8-x^3*(z+x)^5;
197  ideal adj = adjointIdeal(f);
198  def Rn = invertBirMap(adj,ideal(f));
199  setring(Rn);
200  J;
201  psi;
202}
203
204///////////////////////////////////////////////////////////////////////////////
205static proc checkAssumptions(poly f)
206"USAGE:  checkAssumptions(f); f poly
207RETURN:  1 if assumptions are satisfied, 0 otherwise.@*
208         Assumptions checked are: basering is polynomial ring in 3 variables
209         with coefficients in Q, f is homogeneous and absolutely irreducible
210"
211{
212  def Roriginal = basering;
213  list rl = ringlist(Roriginal);
214  rl[3] = list(list("dp",1:3),list("C",0));
215  if(size(rl[1])>1){ERROR("ground field is not Q");}
216  if(rl[1]!=0){ERROR("ground field is not Q");}
217  if(nvars(Roriginal)!=3)
218  { ERROR("not a projective plane curve: wrong number of variables"); }
219  if(homog(f)==0)
220  { ERROR("not a projective plane curve: polynomial is not homogeneous"); }
221  def RP2 = ring(rl);
222  setring RP2;
223  poly f = fetch(Roriginal,f);
224  if(isIrreducible(f)==0){ERROR("curve is not absolutely irreducible");}
225  setring Roriginal;
226}
227
228///////////////////////////////////////////////////////////////////////////////
229proc paraPlaneCurve(poly f, list #)
230"USAGE:  paraPlaneCurve(f [, s]); f poly , s optional string@*
231         optional string s can be: @*
232         'normal': compute integral basis via normalization. @*
233         'local':  make local analysis of singularities first and apply
234                   normalization separately. @*
235         The default is 2. @*
236ASSUME:  The basering must be a polynomial ring in three variables, say x,y,z,
237         with coefficients in Q. @*
238         The polynomial f must be homogeneous and absolutely irreducible. @*
239         The curve C = {f = 0} must be rational, i.e., have geometric genus 0
240         (see @ref{genus}). @*
241         These conditions will be checked automatically.
242RETURN:  ring with an ideal PARA which contains a rational parametrization of
243         the rational plane curve given by f; the ground field of the returned
244         polynomial ring is either Q or some algebraic extension Q(a); PARA
245         consists of three generators that parametrize the three coordinates
246         of the rational curve
247THEORY:  After a first step, realized by a projective automorphism in the
248         procedure adjointIdeal, C satisfies: @*
249- C does not have singularities at infinity z=0. @*
250- C does not contain the point (0:1:0) (that is, the dehomogenization of f
251  with respect to z is monic as a polynomial in y). @*
252Considering C in the chart z<>0, the algorithm regards x as transcendental
253and y as algebraic and computes an integral basis in C(x)[y] of the integral
254closure of C[x] in C(x,y) using the normalization algorithm from @ref{normal_lib}:
255see @ref{integralbasis_lib}. In a future edition of the library, also van Hoeij's
256algorithm for computing the integral basis will be available. @*
257From the integral basis, the adjoint ideal is obtained by linear algebra.
258Alternatively, the algorithm starts with a local analysis of the singular
259locus of C. Then, for each  primary component of the singular locus which
260does not correspond to ordinary multiple points or cusps, the integral
261basis algorithm is applied separately. The ordinary multiple points and
262cusps, in turn, are addressed by a straightforward direct algorithm. The
263adjoint ideal is obtained by intersecting all ideals obtained locally.
264The local variant of the algorithm is used by default. @*
265The linear system corresponding to the adjoint ideal maps the curve
266birationally to a rational normal curve in P^(n-2). @*
267Iterating the anticanonical map, the algorithm projects the rational normal
268curve to PP1 for n odd resp. to a conic C2 in PP2 for n even. @*
269In case n is even, the algorithm tests whether there is a rational point on C2
270and if so gives a parametrization of C2 which is defined over Q. Otherwise the
271parametrization is defined over a quadratic field extension of Q. @*
272By inverting the birational map of C to PP1 resp. to C2, a parametrization of
273C is obtained (defined over Q or the quadratic field extension).
274KEYWORDS: rational curves, rational parametrization of rational curves.
275EXAMPLE:  example paraPlaneCurve; shows an example
276"
277{
278  int choice = 2;
279  if (size(#) != 0)
280  {
281    if (typeof(#[1]) == "string")
282    {
283      string s = string(#[1]);
284      if (s == "normal") { choice = 1; }
285      else
286      {
287        if (s == "local") { choice = 2; }
288        else { ERROR("expected optional argument to be either"
289                   + " 'local' or 'normal'"); }
290      }
291    }
292    else { ERROR("expected optional argument to be a string"); }
293  }
294  def Roriginal = basering;
295  /*checking assumptions and handling the conic case*/
296  checkAssumptions(f);
297  list rl = ringlist(Roriginal);
298  rl[2] = list("x","y","z");
299  rl[3] = list(list("dp",1:3),list("C",0));
300  def RP2 = ring(rl);
301  setring RP2;
302  poly f = fetch(Roriginal,f);
303  int d = deg(f);
304  if(d==2)
305    {
306      def RP1 = paraConic(f); // (ring, PARACONIC)
307      setring RP1;
308      ideal PARA = PARACONIC;
309      export(PARA);
310      "// 'paraPlaneCurve' created a ring together with an ideal PARA.";
311      "// Supposing you typed, say,  def RP1 = paraPlaneCurve(f);";
312      "// you may access the ideal by typing";
313      "//      setring RP1; PARA;";
314      return(RP1);
315    }
316  int k;
317  /*the adjoint ideal*/
318  ideal AI = adjointIdeal(f,list(choice,"rattestyes/firstchecksdone"));
319  /*rattestyes -> causes error message if curve is not rational*/
320  /*firstchecksdone -> prevents that check of assumptions will be done again*/
321  /*mapping the curve to a rational normal curve V(RNC) in Proj(Rrnc) via AI*/
322  def Rrnc = mapToRatNormCurve(f, AI);  // ring containing ideal RNC
323  setring Rrnc;
324  int m = d-1;  // the size of AI
325  /*the odd dimensional case*/
326  if((d mod 2) == 1)
327    {
328      /*mapping the rational normal curve to P^1 creating PHI*/
329      ideal PHI = rncItProjOdd(RNC);
330      /*composing the maps AI and PHI*/
331      def Rbig = Rrnc + RP2;
332      setring Rbig;
333      ideal AI = imap(RP2,AI);
334      ideal PROJ = imap(Rrnc,PHI);
335      for(k = 1; k <= m; k++)
336        {
337          PROJ = subst(PROJ,var(k),AI[k]);
338        }
339      setring RP2;
340      ideal If = ideal(fetch(Roriginal,f));
341      ideal PROJ = imap(Rbig,PROJ);
342      /*inverting the composed map to psi*/
343      def rp1 = invertBirMap(PROJ,If); // ring containing ideal psi
344      setring rp1;
345      list rl1 = ringlist(rp1);
346      rl1[2] = list("s","t");
347      def RP1 = ring(rl1);
348      setring RP1;
349      ideal PARA = fetch(rp1,psi);
350      export(PARA);
351      "// 'paraPlaneCurve' created a ring together with an ideal PARA.";
352      "// Supposing you typed, say,  def RP1 = paraPlaneCurve(f);";
353      "// you may access the ideal by typing";
354      "//      setring RP1; PARA;";
355      return(RP1);
356    }
357  /*the even dimensional case*/
358  /*mapping the rational normal curve to a CONIC in P^2* creating PHI*/
359  def RP2conic = rncItProjEven(RNC);  // exports PHI, returns ring
360                                      // containing CONIC
361  setring RP2conic;
362  /*mapping the conic to P^1 via pencil defined by Ipoint*/
363  def RP2conicNew = projConic(CONIC);  // ring containing ideal Ipoint
364                               // possibly defined over algebraic number field
365                               // variables u,v,w
366  /*composing the maps AI and PHI and Ipoint in two steps*/
367  def Rbig = RP2conicNew + Rrnc;
368  setring Rbig;
369  ideal PHI = imap(Rrnc,PHI);
370  ideal PROJ = imap(RP2conicNew,Ipoint);
371  for(k = 1; k <= 3; k++)
372    {
373      PROJ = subst(PROJ,var(k),PHI[k]);
374    }
375  ideal AI = fetch(RP2,AI);
376  for(k = 1; k <= m; k++)
377    {
378       PROJ = subst(PROJ,var(k+3),AI[k]);
379    }
380  setring RP2conicNew;
381  ideal If = ideal(fetch(Roriginal,f));
382  ideal PROJ = imap(Rbig,PROJ);
383  /*inverting the composed map to psi*/
384  def rp1 = invertBirMap(PROJ,If); // (ring, (J,psi))
385  setring rp1;
386  list rl1 = ringlist(rp1);
387  rl1[2] = list("s","t");
388  def RP1 = ring(rl1);
389  setring RP1;
390  ideal PARA = fetch(rp1,psi);
391  export(PARA);
392  "// 'paraPlaneCurve' created a ring together with an ideal PARA.";
393  "// Supposing you typed, say,  def RP1 = paraPlaneCurve(f);";
394  "// you may access the ideal by typing";
395  "//      setring RP1; PARA;";
396  return(RP1);
397}
398
399example
400{ "EXAMPLE:"; echo=2;
401  ring R = 0,(x,y,z),dp;
402  poly f1 = 1/2*x^5+x^2*y*z^2+x^3*y*z+1/2*x*y^2*z^2-2*x*y^3*z+y^5;
403  def Rp1 = paraPlaneCurve(f1);
404  setring Rp1;
405  PARA;
406  setring R;
407  poly f2 = x6+3x4y2+3x2y4+y6-4x4z2-34x3yz2-7x2y2z2+12xy3z2+6y4z2;
408  f2 = f2+36x2z4+36xyz4+9y2z4;
409  def Rp2 = paraPlaneCurve(f2);
410  setring Rp2;
411  PARA;
412}
413
414///////////////////////////////////////////////////////////////////////////////
415static proc findCoordChange(poly f, ideal JAC)
416"USAGE:  findCoordChange(f, JAC); f poly, JAC ideal.
417ASSUME:  The polynomial f is homogeneous in three variables, JAC is
418         the Jacobi ideal of f.
419RETURN:  intvec, say a,b,c. After the coordinate change
420         var(3) --> a*var(1)+b*var(2)+c*var(3), the curve {f=0}
421         has no singularities at infinity {var(3)=0}.
422"
423{
424  int h = 2;
425  int a,b,c;
426  ideal Jfinfty;
427  while(h)
428    {
429      c = 1;
430      while(c<=h)
431         {
432           b = 0;
433           while(b<=(h-c))
434              {
435                a = h-b-c;
436                Jfinfty = JAC,a*var(1)+b*var(2)+c*var(3);
437                if(dim(std(Jfinfty)) == 0)
438                  {
439                    return(a,b,c);
440                  }
441                b = b+1;
442              }
443            c = c+1;
444          }
445       h = h+1;
446    }
447}
448
449///////////////////////////////////////////////////////////////////////////////
450proc adjointIdeal(poly f, list #)
451"USAGE:  adjointIdeal(f [, choices]); f polynomial in three variables, choices
452         optional list consisting of one integer or of one string or of one
453         integer followed by one string. @*
454         Optional integer can be: @*
455         1: compute integral basis via normalization. @*
456         2: make local analysis of singularities first and apply normalization
457            separately. @*
458         3: normalization via ideal quotient. @*
459         4: normalization via local ideal quotient. @*
460         The default is 2. @*
461         Optional string may contain substrings: @*
462         - rattestyes -> causes error message if curve is not rational. @*
463         - firstchecksdone -> prevents that check of assumptions will be done
464           more than once.
465ASSUME:  The basering must be a polynomial ring in three variables, say x,y,z,
466         with coefficients in Q. @*
467         The polynomial f must be homogeneous and absolutely irreducible.@*
468         All these conditions will be checked automatically.@*
469RETURN:  ideal, the adjoint ideal of the curve defined by f.
470THEORY:  Considering C in the chart z<>0, the algorithm regards x as
471transcendental and y as algebraic and computes an integral basis in C(x)[y] of
472the integral closure of C[x] in C(x,y) using the normalization algorithm
473from @ref{normal_lib}: see @ref{integralbasis_lib}. In a future edition of the library,
474also van Hoeij's algorithm for computing the integral basis will be available.@*
475From the integral basis, the adjoint ideal is obtained by linear algebra.
476Alternatively, the algorithm starts with a local analysis of the singular
477locus of C. Then, for each  primary component of the singular locus which
478does not correspond to ordinary multiple points or cusps, the integral
479basis algorithm is applied separately. The ordinary multiple points and
480cusps, in turn, are addressed by a straightforward direct algorithm. The
481adjoint ideal is obtained by intersecting all ideals obtained locally.
482The local variant of the algorithm is used by default. @*
483KEYWORDS: integral basis; normalization.
484EXAMPLE: example adjointIdeal; shows an example
485"
486{
487  list choices = #;
488  if(size(#)==0)
489    {
490      checkAssumptions(f);
491      choices = list(2, "rattestno");
492    }
493  if(size(#)==1)
494    {
495      if(typeof(choices[1])=="int")
496        {
497          checkAssumptions(f);
498          choices = list(choices[1], "rattestno");
499        }
500      else
501         {
502           if(not(find(choices[1], "firstchecksdone")))
503             {
504               checkAssumptions(f);
505             }
506           else
507             {
508               choices = list(2, choices[1]);
509             }
510        }
511    }
512  if(size(#) == 2)
513    {
514      if(not(find(choices[2],"firstchecksdone")))
515        {
516          checkAssumptions(f);
517        }
518    }
519  ideal JAC = diff(maxideal(1),ideal(f));
520  ideal Jfinfty = JAC,var(3);
521  /*applying a change of coordinates if (f=0) has singularities at infinity*/
522  int bb1;
523  if(dim(std(Jfinfty)) >= 1)
524    {
525       bb1 = 1;
526       int a,b,c  = findCoordChange(f,JAC);
527       f = subst(f,var(3),var(3)-a/c*var(1)-b/c*var(2));
528    }
529  /*applying a change of coordinates if the point (0:1:0) lies on the curve*/
530  matrix co = coeffs(f,var(2));
531  int bb2 = ((size(co)-1) != deg(f));
532  if(bb2)
533    {
534       co = coeffs(f,var(1));
535       int bb2x = ((size(co)-1) == deg(f));
536       if(bb2x)
537         {
538           map perm = basering, var(2), var(1), var(3);
539           f = perm(f);
540         }
541       else
542         {
543           f = subst(f,var(1),var(1)+var(2));
544         }
545    }
546  co = coeffs(f,var(2));
547  f = f/co[size(co),1];
548  /*the actual computation*/
549  ideal AI = adjointIdealAtWork(f,choices);
550  /*reversing the changes of coordinates if needed*/
551  if(bb2)
552    {
553       if(bb2x)
554         {
555           map perm = basering, var(2), var(1), var(3);
556           AI = mstd(perm(AI))[2];
557         }
558       else
559         {
560           AI = mstd(substitute(AI,var(1),var(1)-var(2)))[2];
561         }
562
563    }
564  if(bb1==1)
565   {
566      AI = mstd(substitute(AI,var(3),var(3)+a/c*var(1)+b/c*var(2)))[2];
567    }
568  return(AI);
569}
570
571example
572{ "EXAMPLE:"; echo=2;
573  ring R = 0,(x,y,z),dp;
574  poly f = y^8-x^3*(z+x)^5;
575  adjointIdeal(f);
576}
577
578///////////////////////////////////////////////////////////////////////////////
579static proc adjointIdealAtWork(poly f, list choices)
580"USAGE:  adjointIdealAtWork(f, choices); f polynomial in three variables,
581         choices list consisting of one integer followed by one string. @*
582         integer can be: @*
583         1: compute integral basis via normalization. @*
584         2: make local analysis of singularities first and apply normalization
585            separately. @*
586         3: normalization via ideal quotient. @*
587         4: normalization via local ideal quotient. @*
588         The default is 2. @*
589         string  may contain substring: @*
590         - rattestyes -> causes error message if curve is not rational. @*
591ASSUME:  The basering must be a polynomial ring in three variables, say x,y,z,
592         with coefficients in Q. @*
593         The polynomial f must be homogeneous and absolutely irreducible. @*
594         Its dehomogenization with respect to the third variable must be monic
595         as a polynomial in the second variable (that is, C does not contain
596         the point (0:1:0)).@*
597         The curve C is not allowed to have singularities
598         at infinity (z = 0). @*
599RETURN:  ideal, the adjoint ideal of the curve defined by f.
600"
601{
602  def Roriginal = basering;
603  list rl = ringlist(Roriginal);
604  rl[3] = list(list("dp",1:nvars(Roriginal)),list("C",0));
605  def RP2 = ring(rl);
606  setring RP2;
607  poly f = imap(Roriginal,f);
608  poly dhf = subst(f,var(3),1);
609  int n = deg(f);
610  if((choices[1]==1) || (choices[1]==3)) // no local analysis of singularities
611    {
612      ideal AI;
613      if (choices[1]==1)
614      { AI = adjointIdealIB(f,insert(choices,ideal(0),size(choices))); }
615      else
616      { AI = adjointIdealIQ(f,insert(choices,ideal(0),size(choices))); }
617      AI = homog(std(AI),var(3));
618      AI = sat(AI, maxideal(1))[1];
619      AI = minbase(AI);
620      setring Roriginal;
621      return(imap(RP2,AI));
622    }
623  list LL = geomGenusLA(f);  // local analysis of singularities
624  int sizeLL2 = size(LL[2]);
625  int sizeLL3 = size(LL[3]);
626  list LL3 = LL[3];
627  ideal LL4 = LL[4];
628  if((LL[1]!=0) && (find(choices[2],"rattestyes")))
629  { ERROR("not a rational curve"); }
630  if((LL[2]==0) && (sizeLL3==0)  && (LL4==1))  // smooth case
631    {
632      setring Roriginal;
633      return(ideal(1));
634    }
635  int j,k;
636  list rl = ringlist(RP2);
637  rl[2] = list(var(1), var(2));
638  rl[3] = list(list("dp",1:2),list("C",0));
639  def Rdummy = ring(rl);
640  ideal B;
641  if(sizeLL3==0){B = 1;} // no ordinary multiple points
642                               // (except possibly nodes)
643  else                         // there are ordinary multiple points
644                               // (other than nodes)
645    {
646      setring Rdummy;
647      list OMP = imap(RP2,LL3);
648      int ub;
649      for(k=1;k<=size(OMP);k++)
650        {
651          if(OMP[k][1]>ub)
652            {
653              ub = OMP[k][1];
654            }
655        }
656     int lb = ub;
657     for(k=1;k<=size(OMP);k++)
658        {
659          if(OMP[k][1]<lb)
660            {
661              lb = OMP[k][1];
662            }
663        }
664     for(k=lb;k<=ub;k++)
665       {
666         ideal A(k) = 1;
667       }
668     for(k=1;k<=size(OMP);k++)
669       {
670         A(OMP[k][1]) = intersect(A(OMP[k][1]), OMP[k][2]);
671       }
672     int i = ub;
673     setring RP2;
674     for(k=lb;k<=ub;k++)
675       {
676         ideal A(k) = homog(std(fetch(Rdummy,A(k))),var(3));
677       }
678     B = maxideal(n-i);
679     ideal A;
680     while(i>=lb)
681        {
682          A = A(i)**(i-1);
683          j=1;
684          while(j<=ncols(A))
685            {
686              if(deg(A[j]>(n-2)))
687                {
688                  A = sat(A, maxideal(1))[1];
689                  break;
690                 }
691              j = j+1;
692            }
693          B = intersect(B,A);
694          i = i-1;
695       }
696    }  //end else
697  B = intersect(B,homog(std(LL4),var(3)));  // add nodes and cusps
698  if(sizeLL2==0)  // ordinary multiple points plus cusps only
699    {
700      ideal AI = sat(B, maxideal(1))[1];
701      AI = minbase(AI);
702      setring Roriginal;
703      return(imap(RP2,AI));
704    }
705  setring Rdummy;
706  poly f = imap(RP2,dhf);
707  ideal SL = jacob(f),f;
708  SL = sat(SL, fetch(RP2,LL4))[1];
709  if(sizeLL3!=0)
710    {
711      for(k=lb;k<=ub;k++)
712        {
713          SL = sat(SL, A(k))[1];
714        }
715    }
716  list PD = primdecGTZ(SL);  // could be made faster -- see minAssGTZ
717                             // in deltaLocMod -- only one PD needed
718  int pd = size(PD);
719  setring RP2;
720  list PD = imap(Rdummy,PD);
721  ideal AI = 1;
722  for(k=1;k<=pd;k++)
723    {
724      if (choices[1]==2)
725      { AI = intersect(AI,adjointIdealIB(f,insert(choices,PD[k][1],
726                                                 size(choices)))); }
727      else
728      { AI = intersect(AI,adjointIdealIQ(f,insert(choices,PD[k][1],
729                                                 size(choices)))); }
730    }
731  AI = homog(std(AI),var(3));
732  AI = intersect(AI,B);
733  AI = sat(AI, maxideal(1))[1];
734  AI = minbase(AI);
735  setring Roriginal;
736  return(imap(RP2,AI));
737}
738
739///////////////////////////////////////////////////////////////////////////////
740static proc adjointIdealIB(poly f, list choices)
741"USAGE:  adjointIdealIB(f, choices); f polynomial in three variables, choices
742         list consisting of one integer followed by one string followed by one
743         ideal. @*
744         integer can be: @*
745         1, 2 : compute integral basis via normalization @*
746         The default is 2. @*
747         string  may contain substring: @*
748         - rattestyes -> causes error message if curve is not rational. @*
749         ideal serves as input for @ref integralBasis.
750ASSUME:  The basering must be a polynomial ring in three variables, say x,y,z,
751         with coefficients in Q. @*
752         The polynomial f must be homogeneous and absolutely irreducible.@*
753         Its dehomogenization with respect to the third variable must be monic
754         as a polynomial in the second variable (that is, C does not contain
755         the point (0:1:0)).@*
756         The curve C is not allowed to have singularities
757         at infinity (z = 0). @*
758RETURN:  ideal containing the adjoint ideal of the curve defined by f. @*
759"
760{
761"#############################################";
762int totalTime = timer;
763  poly dhf = subst(f,var(3),1);
764  def Roriginal = basering;
765  list rl = ringlist(Roriginal);
766  rl[2] = list(var(1), var(2));
767  rl[3] = list(list("dp",1:2),list("C",0));
768  def Rdummy = ring(rl);
769  setring Rdummy;
770  poly f = imap(Roriginal,dhf);
771  poly d2f = diff(f,var(2));
772  list DATA = imap(Roriginal,choices);
773  /* Creating rings for later use */
774  list rl = ringlist(Rdummy);
775  rl[2] = list(var(2), var(1));
776  rl[3] = list(list("lp",1:2),list("C",0));
777  def Rred = ring(rl);   // make var(2) > var(1)
778  rl = ringlist(Rdummy);
779  rl[1] = list(0,list(var(1)),list(list("dp",1)),ideal(0));
780  rl[2] = list(var(2));
781  rl[3] = list(list("dp",1),list("C",0));
782  def QF = ring(rl);   // make var(1) transcendental
783  list LIntB;
784int t = timer;
785  if(DATA[1] <= 4)  // use normalization algorithm
786    {
787      LIntB = integralBasis(f, 2, list(list("inputC", DATA[3]),"isIrred"));
788    }
789  else                                 // use van Hoeij's algorithm
790    {
791      LIntB = integralBasisVH(f,DATA[3],2);  // van Hoeij in future version
792                                             // used when DATA[1] = 5
793    }
794  if(find(DATA[2],"rattestyes") && (DATA[3]==0))
795    {
796      setring Roriginal;
797      int gg = geomGenusIB(f,imap(Rdummy, LIntB));
798      if(gg!=0){ERROR("not a rational curve");}
799      setring Rdummy;
800    }
801t = timer - t;
802"### adjointIdealIB, IB-time:", t, "msec";
803  int i,j,k,l;
804  ideal IB = LIntB[1];
805  poly d = LIntB[2];
806  int sL=size(IB);
807  setring Rred;
808  ideal IB = imap(Rdummy,IB);
809  ideal fred = std(imap(Rdummy,f));
810  IB = reduce(IB,fred);
811  matrix M = coeffs(IB,var(1));
812  setring QF;
813  matrix M = imap(Rred,M);
814  poly d = imap(Rdummy,d);
815  M=1/d*M;
816t = timer;
817  list LUM = ludecomp(M);
818t = timer - t;
819"### adjointIdealIB, time for 1st ludecomp:", t, "msec";
820  list LS;
821  matrix dummyvector[sL][1];
822  matrix Gij[sL][sL];
823  matrix Tr[sL][sL];
824  setring Rred;
825  poly eiej;
826  list Iij, empty;
827  matrix Gij[sL][sL];
828t = timer;
829  for(i = 1; i <= sL; i++)
830     {
831       for(j = i; j <= sL; j++)
832          {
833            setring Rred;
834            Gij = 0;
835            eiej = IB[i]*IB[j];
836            Iij=empty;
837            for(k = 1; k <= sL; k++)
838               {
839                  Iij[k] = reduce(eiej*IB[k],fred);
840               }
841            Gij = coeffs(ideal(Iij[1..sL]),var(1));
842            setring QF;
843            Gij = imap (Rred, Gij);
844            for(k = 1; k <= sL; k++)
845               {
846                  dummyvector = Gij[1..sL,k];
847                  LS = lusolve(LUM[1], LUM[2], LUM[3], dummyvector);
848                  Tr[i,j] = Tr[i,j] + 1/d^3*LS[2][k,1];
849               }
850          }
851     }
852  for(i = 1; i <= sL; i++)
853     {
854       for(j = 1; j < i; j++)
855          {
856             Tr[i,j] = Tr[j,i];
857          }
858     }
859t = timer - t;
860"### adjointIdealIB, time for linear solving:", t, "msec";
861t = timer;
862  LUM = ludecomp(Tr);
863t = timer - t;
864"### adjointIdealIB, time for 2nd ludecomp:", t, "msec";
865  setring Rred;
866  poly d2f = imap(Rdummy,d2f);
867  IB = d2f*IB;
868  IB = reduce(IB,fred);
869  setring QF;
870  matrix IB = transpose(matrix(imap(Rred,IB)));
871  IB = 1/d*IB;
872  LS = lusolve(LUM[1], LUM[2], LUM[3], IB);
873  ideal LL = ideal(LS[2]);
874  setring Roriginal;
875  ideal AI = imap(QF,LL);
876totalTime = timer - totalTime;
877"### adjointIdealIB, total time:", totalTime, "msec";
878  return(AI);
879}
880
881///////////////////////////////////////////////////////////////////////////////
882static proc adjointIdealIQ(poly f, list choices)
883{
884"#############################################";
885int totalTime = timer;
886  poly dhf = subst(f,var(3),1);
887  def Roriginal = basering;
888  list rl = ringlist(Roriginal);
889  rl[2] = list(var(1), var(2));
890  rl[3] = list(list("dp",1:2),list("C",0));
891  def Rdummy = ring(rl);
892  setring Rdummy;
893  list DATA = imap(Roriginal,choices);
894  poly f = imap(Roriginal,dhf);
895  list LIntB;
896int t = timer;
897  if(DATA[1] <= 4)  // use normalization algorithm
898    {
899      LIntB = integralBasis(f, 2, list(list("inputC", DATA[3]),"isIrred"));
900    }
901  else                                 // use van Hoeij's algorithm
902    {
903      LIntB = integralBasisVH(f,DATA[3],2);  // van Hoeij in future version
904                                             // used when DATA[1] = 5
905    }
906  if(find(DATA[2],"rattestyes") && (DATA[3]==0))
907    {
908      setring Roriginal;
909      int gg = geomGenusIB(f,imap(Rdummy, LIntB));
910      if(gg!=0){ERROR("not a rational curve");}
911      setring Rdummy;
912    }
913t = timer - t;
914"### adjointIdealIQ, IB-time:", t, "msec";
915  int i,j,k,l;
916  ideal IB = LIntB[1];
917  poly d = LIntB[2];
918  ideal fd = f, d;
919  ideal IBf = IB, f;
920t = timer;
921  ideal AI = quotient(fd, IBf);
922t = timer - t;
923"### adjointIdealIQ, quot-time:", t, "msec";
924//"#### IB:"; IB;
925//"#### d:", d;
926  setring Roriginal;
927  ideal AI = imap(Rdummy,AI);
928//"#### AI:"; AI;
929totalTime = timer - totalTime;
930"### adjointIdealIQ, total time:", totalTime, "msec";
931  return(AI);
932}
933
934///////////////////////////////////////////////////////////////////////////////
935proc mapToRatNormCurve(poly f, ideal AI)
936"USAGE:  mapToRatNormCurve(f, AI); f polynomial, AI ideal
937ASSUME:  The polynomial f is homogeneous in three variables and absolutely
938         irreducible.
939         The plane curve C defined by f is rational.
940         The ideal AI is the adjoint ideal of C.
941RETURN:  ring with an ideal RNC.
942EXAMPLE: example mapToRatNormCurve; shows an example
943"
944{
945   int n = size(AI);
946   int k;
947   //if(n!=deg(f)-1){ERROR("not a rational curve");}
948   def Roriginal = basering;
949   ideal IC = f;
950   list rl = ringlist(Roriginal);
951   /* begin workaround elimination*/
952   for(k = 1; k <= 3; k++)
953     {
954        rl[2][k] = "x("+string(k)+")";
955     }
956   for(k = 1; k <= n; k++)
957     {
958        rl[2][k+3] = "y("+string(k)+")";
959     }
960   rl[3]= list(list("dp",1:(3+n)),list("C",0));
961   def Relim = ring(rl);
962   setring Relim;
963   ideal IC = fetch(Roriginal,IC);
964   ideal AI = fetch(Roriginal,AI);
965   ideal J;
966   J = IC;
967   for(k=1;k<=n;k++)
968     {
969       J=J,var(k+3)-AI[k];
970     }
971   ideal SJ = std(J);
972   intvec HJ = hilb(SJ,1);
973   ideal RNC = eliminate(J,x(1)*x(2)*x(3),HJ);
974   list rl = ringlist(Relim);
975   list rl2 = rl[2];
976   rl[2] = list(rl2[4..n+3]);
977   rl[3]= list(list("dp",1:n),list("C",0));
978   def Rtarget = ring(rl);
979   setring Rtarget;
980   ideal RNC = imap(Relim,RNC);
981   /* end workaround elimination*/
982   export(RNC);
983   int p = printlevel - voice + 3;
984   dbprint(p,"//'mapToRatNorm' created a ring together with an ideal RNC.");
985   dbprint(p,"// Supposing you typed, say,  def RPn = mapToRatNorm(f,AI);");
986   dbprint(p,"// you may access the ideal by typing");
987   dbprint(p,"//      setring RPn; RNC;");
988   return(Rtarget);
989}
990
991example
992{ "EXAMPLE:"; echo=2;
993  ring R = 0,(x,y,z),dp;
994  poly f = y^8-x^3*(z+x)^5;
995  ideal adj = adjointIdeal(f);
996  def Rn = mapToRatNormCurve(f,adj);
997  setring(Rn);
998  RNC;
999}
1000
1001///////////////////////////////////////////////////////////////////////////////
1002proc rncAntiCanonicalMap(ideal I)
1003"USAGE:  rncAntiCanonicalMap(I); I ideal
1004ASSUME:  I is a homogeneous ideal in the basering
1005         defining a rational normal curve C in PP^n.
1006NOTE:   The procedure will fail or give a wrong output if I is not the
1007        ideal of a rational normal curve.
1008RETURN:  ideal defining the anticanonical map  C --> PP^(n-2). @*
1009         Note that the entries of the ideal should be considered as
1010         representatives of elements in R/I, where R is the basering.
1011THEORY:  The anti-canonical map of a rational normal curve
1012         maps C isomorpically to a rational normal curve in PP^(n-2).
1013KEYWORDS: rational normal curve, projection.
1014EXAMPLE: example rncAntiCanonicalMap; shows an example
1015"
1016{
1017  def Roriginal = basering;
1018  list rl = ringlist(Roriginal);
1019  rl[3] = list(list("dp",1:nvars(Roriginal)),list("C",0));
1020  def RoriginalDP = ring(rl);
1021  setring RoriginalDP;
1022  ideal I = imap(Roriginal,I);
1023  int cc = nvars(RoriginalDP)-2;
1024  module AKD = Ext_R(cc,I);
1025  qring qI = std(I);
1026  matrix AKD = imap(RoriginalDP,AKD);
1027  AKD = syz(transpose(AKD));
1028  ideal PR = submat(AKD,1..nrows(AKD),1);
1029  setring Roriginal;
1030  return(imap(qI,PR));
1031}
1032
1033example
1034{ "EXAMPLE:"; echo=2;
1035  ring R = 0,(x,y,z),dp;
1036  poly f = y^8-x^3*(z+x)^5;
1037  ideal adj = adjointIdeal(f);
1038  def Rn = mapToRatNormCurve(f,adj);
1039  setring(Rn);
1040  RNC;
1041  rncAntiCanonicalMap(RNC);
1042}
1043
1044
1045///////////////////////////////////////////////////////////////////////////////
1046proc rncItProjOdd(ideal I)
1047"USAGE:  rncItProjOdd(I); I ideal
1048ASSUME:  I is a homogeneous ideal in the basering with n+1 variables
1049         defining a rational normal curve C in PP^n with n odd.
1050NOTE:    The procedure will fail or give a wrong output if I is not the
1051         ideal of a rational normal curve. It will test whether n is odd.
1052RETURN:  ideal PHI defining an isomorphic projection of C to PP^1.@*
1053         Note that the entries of PHI should be considered as
1054         representatives of elements in R/I, where R is the basering.
1055THEORY:  We iterate the procedure @ref{rncAntiCanonicalMap} to obtain PHI.
1056KEYWORDS: rational normal curve, projection.
1057SEE ALSO: rncItProjEven.
1058EXAMPLE: example rncItProjOdd; shows an example
1059"
1060{
1061  int n = nvars(basering);
1062  if((n mod 2) == 1){ERROR("Pn has even dimension");}
1063  def Roriginal = basering;
1064  list rlo = ringlist(Roriginal);
1065  rlo[3]= list(list("dp",1:n),list("C",0));
1066  int k;
1067  for(k = 1; k <= n; k++)
1068    {
1069      rlo[2][k] = "z("+string(k)+")";
1070    }
1071  def RoriginalCopy = ring(rlo);
1072  for(k = 1; k <= n; k++)
1073    {
1074      rlo[2][k] = "y("+string(k)+")";
1075    }
1076  def Rold = ring(rlo);
1077  setring RoriginalCopy;
1078  ideal PHI  = maxideal(1);
1079  setring Rold;
1080  ideal J = fetch(Roriginal,I);
1081  list rl2;
1082  def Rnew;
1083  def Rbig;
1084  def Relim;
1085  intvec HJJ;
1086  while(n>2)
1087     {
1088        ideal PR = rncAntiCanonicalMap(J);
1089        list rl = ringlist(Rold);
1090        Rbig = Rold + RoriginalCopy;
1091        setring Rbig;
1092        ideal PHI = imap(RoriginalCopy,PHI);
1093        ideal dummy = imap(Rold,PR);
1094        for(k = 1; k <= n; k++)
1095          {
1096             dummy = subst(dummy,var(k),PHI[k]);
1097          }
1098        setring RoriginalCopy;
1099        PHI = imap(Rbig,dummy);
1100        /* begin workaround elimination*/
1101        setring Rold;
1102        for(k = 1; k <= n; k++)
1103          {
1104            rl[2][k] = "x("+string(k)+")";
1105          }
1106        for(k = 1; k <= n-2; k++)
1107          {
1108            rl[2][k+n] = "y("+string(k)+")";
1109          }
1110        rl[3]= list(list("dp",1:(2*n-2)),list("C",0));
1111        Relim = ring(rl);
1112        setring Relim;
1113        ideal J = fetch(Rold,J);
1114        ideal PR = fetch(Rold,PR);
1115        ideal JJ = J;
1116        poly pvar=1;
1117        for(k = 1; k <= n; k++)
1118          {
1119            pvar = pvar*var(k);
1120          }
1121        for(k=1;k<=n-2;k++)
1122          {
1123            JJ=JJ,var(k+n)-PR[k];
1124          }
1125        ideal SJJ = std(JJ);
1126        HJJ = hilb(SJJ,1);
1127        J = eliminate(JJ,pvar,HJJ);
1128        list rl = ringlist(Relim);
1129        rl2 = rl[2];
1130        rl[2] = list(rl2[n+1..2*n-2]);
1131        rl[3]= list(list("dp",1:(n-2)),list("C",0));
1132        Rnew = ring(rl);
1133        setring Rnew;
1134        ideal J = imap(Relim,J);
1135        /* end workaround elimination*/
1136        Rold = Rnew;
1137        setring Rold;
1138        n = n-2;
1139     }
1140  setring Roriginal;
1141  return(fetch(RoriginalCopy,PHI));
1142}
1143
1144example
1145{ "EXAMPLE:"; echo=2;
1146  ring R = 0,(x,y,z),dp;
1147  poly f = -x7-10x5y2-10x4y3-3x3y4+8x2y5+7xy6+11y7+3x6+10x5y +30x4y2
1148           +26x3y3-13x2y4-29xy5-33y6-3x5-20x4y-33x3y2-8x2y3+37xy4+33y5
1149           +x4+10x3y+13x2y2-15xy3-11y4;
1150  f = homog(f,z);
1151  ideal adj = adjointIdeal(f);
1152  def Rn = mapToRatNormCurve(f,adj);
1153  setring(Rn);
1154  RNC;
1155  rncItProjOdd(RNC);
1156}
1157
1158
1159///////////////////////////////////////////////////////////////////////////////
1160proc rncItProjEven(ideal I)
1161"USAGE:  rncItProjEven(I); I ideal
1162ASSUME:  I is a homogeneous ideal in the basering with n+1 variables
1163         defining a rational normal curve C in PP^n with n even.
1164NOTE:    The procedure will fail or give a wrong output if I is not the
1165         ideal of a rational normal curve. It will test whether n is odd.
1166RETURN:  ring with an ideal CONIC defining a conic C2 in PP^2.@*
1167         In addition, an ideal PHI in the basering defining an isomorphic
1168         projection of C to C2 will be exported.@*
1169         Note that the entries of PHI should be considered as
1170         representatives of elements in R/I, where R is the basering.
1171THEORY:  We iterate the procedure @ref{rncAntiCanonicalMap} to obtain PHI.
1172KEYWORDS: rational normal curve, projection.
1173SEE ALSO: rncItProjOdd.
1174EXAMPLE: example rncItProjEven; shows an example
1175"
1176{
1177  int n = nvars(basering);
1178  if((n mod 2) == 0){ERROR("Pn has odd dimension");}
1179  def Roriginal = basering;
1180  list rlo = ringlist(Roriginal);
1181  rlo[3]= list(list("dp",1:n),list("C",0));
1182  int k;
1183  for(k = 1; k <= n; k++)
1184    {
1185      rlo[2][k] = "z("+string(k)+")";
1186    }
1187  def RoriginalCopy = ring(rlo);
1188  for(k = 1; k <= n; k++)
1189    {
1190      rlo[2][k] = "y("+string(k)+")";
1191    }
1192  def Rold = ring(rlo);
1193  setring RoriginalCopy;
1194  ideal PHI  = maxideal(1);
1195  setring Rold;
1196  ideal J = fetch(Roriginal,I);
1197  list rl2;
1198  def Rnew;
1199  def Rbig;
1200  def Relim;
1201  intvec HJJ;
1202  while(n>3)
1203     {
1204        ideal PR = rncAntiCanonicalMap(J);
1205        list rl = ringlist(Rold);
1206        Rbig = Rold + RoriginalCopy;
1207        setring Rbig;
1208        ideal PHI = imap(RoriginalCopy,PHI);
1209        ideal dummy = imap(Rold,PR);
1210        for(k = 1; k <= n; k++)
1211          {
1212             dummy = subst(dummy,var(k),PHI[k]);
1213          }
1214        setring RoriginalCopy;
1215        PHI = imap(Rbig,dummy);
1216        /* begin workaround elimination*/
1217        setring Rold;
1218        for(k = 1; k <= n; k++)
1219          {
1220            rl[2][k] = "x("+string(k)+")";
1221          }
1222        for(k = 1; k <= n-2; k++)
1223          {
1224            rl[2][k+n] = "y("+string(k)+")";
1225          }
1226        rl[3]= list(list("dp",1:(2*n-2)),list("C",0));
1227        Relim = ring(rl);
1228        setring Relim;
1229        ideal J = fetch(Rold,J);
1230        ideal PR = fetch(Rold,PR);
1231        ideal JJ = J;
1232        poly pvar=1;
1233        for(k = 1; k <= n; k++)
1234          {
1235            pvar = pvar*var(k);
1236          }
1237        for(k=1;k<=n-2;k++)
1238          {
1239            JJ=JJ,var(k+n)-PR[k];
1240          }
1241        ideal SJJ = std(JJ);
1242        HJJ = hilb(SJJ,1);
1243        J = eliminate(JJ,pvar,HJJ);
1244        list rl = ringlist(Relim);
1245        rl2 = rl[2];
1246        rl[2] = list(rl2[n+1..2*n-2]);
1247        rl[3]= list(list("dp",1:(n-2)),list("C",0));
1248        Rnew = ring(rl);
1249        setring Rnew;
1250        ideal J = imap(Relim,J);
1251        /* end workaround elimination*/
1252        Rold = Rnew;
1253        setring Rold;
1254        n = n-2;
1255     }
1256  poly CONIC = J[1];
1257  export(CONIC);
1258  setring Roriginal;
1259  ideal PHI = fetch(RoriginalCopy,PHI);
1260  export(PHI);
1261  return(Rold);
1262}
1263
1264example
1265{ "EXAMPLE:"; echo=2;
1266  ring R = 0,(x,y,z),dp;
1267  poly f = y^8-x^3*(z+x)^5;
1268  ideal adj = adjointIdeal(f);
1269  def Rn = mapToRatNormCurve(f,adj);
1270  setring(Rn);
1271  RNC;
1272  def Rc = rncItProjEven(RNC);
1273  PHI;
1274  setring Rc;
1275  CONIC;
1276}
1277
1278///////////////////////////////////////////////////////////////////////////////
1279static proc geomGenusIB(poly f, list #)
1280"USAGE:  geomGenusIB(f [, L]); f poly, L optional list representing the
1281         integral basis as returned by @ref integralBasisJ.
1282ASSUME:  The basering must be a polynomial ring in three variables, say x,y,z,
1283         with coefficients in Q. @*
1284         The polynomial f must be homogeneous and absolutely irreducible.@*
1285         Its dehomogenization with respect to the third variable must be monic
1286         as a polynomial in the second variable (that is, the curve C = {f = 0}
1287         does not contain the point (0:1:0)).@*
1288         The curve C is not allowed to have singularities
1289         at infinity (z = 0). @*
1290NOTE:    The last two conditions can be met by a suitable change of coordinates in PGL(3)
1291         as applied in the procedure @ref adjointIdeal. The other conditions
1292         can be tested using @ref checkAssumptions.@*
1293RETURN:  int, the geometric genus of C.
1294THEORY:  We compute an integral basis of the integral closure of the coordinate
1295         ring of C and from that the geometric genus.@*
1296KEYWORDS: geometric genus, plane curves.
1297SEE ALSO: genus.
1298"
1299{
1300  int bb = size(#);
1301  poly dhf = subst(f,var(3),1);
1302  def Roriginal = basering;
1303  list rl = ringlist(Roriginal);
1304  rl[2] = list(var(1), var(2));
1305  rl[3] = list(list("dp",1:2),list("C",0));
1306  def Rdummy = ring(rl);
1307  setring Rdummy;
1308  poly f = imap(Roriginal,dhf);
1309  list LIntB;
1310  if(bb == 0)
1311    {
1312      LIntB = integralBasis(f,2,"isIrred");
1313    }
1314  else
1315    {
1316      LIntB = imap(Roriginal,#);
1317    }
1318  ideal IB = LIntB[1];
1319  poly d = LIntB[2];
1320  int ud = deg(d);
1321  int sL = size(IB);
1322  int k;
1323  int gg = (sL-1)*(sL-2)/2-sL*ud;
1324  for(k = 1; k <= sL; k++)
1325     {
1326        gg = gg + deg(gcd(d,IB[k]));
1327     }
1328  setring Roriginal;
1329  return(gg);
1330}
1331
1332
1333///////////////////////////////////////////////////////////////////////////////
1334static proc geomGenusLA(poly F)
1335"USAGE:  geomGenusLA(F); F polynomial
1336ASSUME:  The basering must be a polynomial ring in three variables. @*
1337         The polynomial F must be homogeneous.@*
1338RETURN:  list L:
1339@texinfo
1340@table @asis
1341@item @code{L[1]}; int:
1342         the geometric genus p_g = p_a - delta of the projective
1343         curve C defined by F, where p_a is the arithmetic genus.
1344@item @code{L[2]}; int:
1345         is positive if C has singularities other
1346         than ordinary multiple points.@*
1347@item @code{L[3]}; list:
1348         consists of one list for each primary component
1349         of the singular locus of C which correponds to a set of conjugated
1350         ordinary multiple points. Each list consists of an int, the
1351         multiplicity of the points, and an ideal, the primary component.
1352@end table
1353@end texinfo
1354NOTE:    delta is the sum of all local delta-invariants of the singularities,
1355         i.e. dim(R'/R), R' the normalization of the local ring R of the
1356         singularity. @*
1357SEE ALSO: genus
1358"
1359{
1360   int w = printlevel-voice+2;  // w=printlevel (default: w=0)
1361   int d = deg(F);
1362   def R = basering;
1363   execute("ring S=("+charstr(R)+"),(x,y,t),dp;");
1364   execute("ring C=("+charstr(R)+"),(x,y),ds;");
1365   int genus=(d-1)*(d-2)/2;
1366   if(w>=1){"the arithmetic genus of the plane curve:";genus;pause();}
1367
1368   int delt,deltaloc,deltainf,tau,tauinf,cusps,iloc,iglob,l,nsing,
1369       tauloc,tausing,k,rat,nbranchinf,nbranch,nodes,cuspsinf,nodesinf;
1370   list inv;
1371   execute("ring newR=("+charstr(R)+"),(x,y),dp;");
1372   //the singularities at the affine part
1373   map sigma=R,var(1),var(2),1;
1374   ideal I=sigma(F);
1375
1376   list OMPButNodes;
1377   int sizeOMPButNodes;
1378   int NotOnlyOMPPlusCusps;
1379
1380   ideal I1=jacob(I);
1381   matrix Hess[2][2]=jacob(I1);
1382   ideal ID=I+I1+ideal(det(Hess));//singular locus of I+I1
1383   ideal radID=std(radical(ID));//the non-nodal locus
1384   if(w>=1){"the non-nodal locus:";"";radID;pause();"";}
1385   if(deg(radID[1])==0)
1386   {
1387     ideal IDsing=1;
1388   }
1389   else
1390   {
1391     ideal IDsing=minor(jacob(ID),2)+radID;//singular locus of ID
1392   }
1393
1394   iglob=vdim(std(IDsing));
1395
1396   ideal radIDsing = 1;
1397
1398   if(iglob!=0)//computation of the radical of IDsing
1399   {
1400      radIDsing=reduce(IDsing,radID);
1401      if(size(radIDsing)==0)
1402      {
1403         radIDsing=radID;
1404         attrib(radIDsing,"isSB",1);
1405      }
1406      else
1407      {
1408         radIDsing=std(radical(IDsing));
1409      }
1410      iglob=vdim(radIDsing);
1411      if((w>=1)&&(iglob))
1412          {"the non-nodal-cuspidal locus:";radIDsing;pause();"";}
1413   }
1414   cusps=vdim(radID)-iglob;
1415
1416   ideal NodesPlusCusps  = radical(sat(I+I1, radIDsing)[1]);
1417
1418   nsing=nsing+cusps;
1419
1420   if(iglob==0)
1421   {
1422      if(w>=1){"             there are only cusps and nodes";"";}
1423      tau=vdim(std(I+jacob(I)));
1424      tauinf=tauinf+tau;
1425      nodes=tau-2*cusps;
1426      delt=nodes+cusps;
1427      nbranch=2*tau-3*cusps;
1428      nsing=nsing+nodes;
1429   }
1430   else
1431   {
1432       if(w>=1){"the non-nodal-cuspidal singularities";"";}
1433       setring C;
1434       ideal I1=imap(newR,radIDsing);
1435       iloc=vdim(std(I1));
1436       if(iglob==iloc)
1437       {
1438          if(w>=1){"only cusps and nodes outside (0,0,1)";}
1439          setring newR;
1440          tau=vdim(std(I+jacob(I)));
1441          tauinf=tauinf+tau;
1442          inv=deltaLocMod(I[1],maxideal(1));
1443          delt=inv[1];
1444          tauloc=inv[2];
1445          nodes=tau-tauloc-2*cusps;
1446          nsing=nsing+nodes;
1447          if(inv[2]!=0)
1448            {
1449              nsing=nsing+1;
1450            }
1451          nbranch=inv[3]+ 2*nodes+cusps;
1452          delt=delt+nodes+cusps;
1453          if((w>=1)&&(inv[2]==0)){"smooth at (0,0,1)";}
1454          if(inv[4]!=0)
1455            {
1456              OMPButNodes = insert(OMPButNodes,list(inv[4],maxideal(1)),
1457                                   sizeOMPButNodes);
1458              sizeOMPButNodes = size(OMPButNodes); // new
1459            }
1460          else
1461            {
1462              NotOnlyOMPPlusCusps = NotOnlyOMPPlusCusps + 1;
1463            }
1464        }
1465        else
1466        {
1467           setring newR;
1468           list pr=minAssGTZ(radIDsing);
1469           if(w>=1){pr;}
1470
1471           for(k=1;k<=size(pr);k++)
1472           {
1473              if(w>=1){nsing=nsing+vdim(std(pr[k]));}
1474              inv=deltaLocMod(I[1],pr[k]);
1475              delt=delt+inv[1];
1476              tausing=tausing+inv[2];
1477              nbranch=nbranch+inv[3];
1478              if(inv[4]!=0)
1479                {
1480                  OMPButNodes = insert(OMPButNodes,list(inv[4],pr[k]),
1481                                       sizeOMPButNodes);
1482                  sizeOMPButNodes = size(OMPButNodes);
1483                }
1484              else
1485                {
1486                  NotOnlyOMPPlusCusps = NotOnlyOMPPlusCusps + 1;
1487                }
1488           }
1489           tau=vdim(std(I+jacob(I)));
1490           tauinf=tauinf+tau;
1491           nodes=tau-tausing-2*cusps;
1492           nsing=nsing+nodes;
1493           delt=delt+nodes+cusps;
1494           nbranch=nbranch+2*nodes+cusps;
1495        }
1496   }
1497   genus=genus-delt-deltainf;
1498   if(w>=1)
1499   {
1500      "The projected plane curve has locally:";"";
1501      "singularities:";nsing;
1502      "branches:";nbranch+nbranchinf;
1503      "nodes:"; nodes+nodesinf;
1504      "cusps:";cusps+cuspsinf;
1505      "Tjurina number:";tauinf;
1506      "Milnor number:";2*(delt+deltainf)-nbranch-nbranchinf+nsing;
1507      "delta of the projected curve:";delt+deltainf;
1508      //"delta of the curve:";p_a-genus;
1509      "genus:";genus;
1510      "====================================================";
1511      "";
1512   }
1513   setring R;
1514   if(sizeOMPButNodes>0)
1515     {
1516       list OMPButNodes = fetch(newR,OMPButNodes);
1517     }
1518   return(list(genus,NotOnlyOMPPlusCusps,OMPButNodes,
1519               fetch(newR,NodesPlusCusps)));
1520}
1521
1522
1523///////////////////////////////////////////////////////////////////////////////
1524static proc deltaLocMod(poly f,ideal singL)
1525"USAGE:  deltaLoc(f,J);  f poly, J ideal
1526ASSUME: f is reduced bivariate polynomial; basering has exactly two variables;
1527        J is irreducible prime component of the singular locus of f (e.g., one
1528        entry of the output of @code{minAssGTZ(I);}, I = <f,jacob(f)>).
1529RETURN:  list L:
1530@texinfo
1531@table @asis
1532@item @code{L[1]}; int:
1533         the sum of (local) delta invariants of f at the (conjugated) singular
1534         points given by J.
1535@item @code{L[2]}; int:
1536         the sum of (local) Tjurina numbers of f at the (conjugated) singular
1537         points given by J.
1538@item @code{L[3]}; int:
1539         the sum of (local) number of branches of f at the (conjugated)
1540         singular points given by J.
1541@item @code{L[3]}; int:
1542         the multiplicity of f at the (conjugated) singular points given by J,
1543         if these are ordinary multiple points, and 0 otherwise.
1544@end table
1545@end texinfo
1546NOTE:    procedure makes use of @code{execute}; increasing printlevel displays
1547         more comments (default: printlevel=0).
1548SEE ALSO: deltaLoc, delta, tjurina
1549KEYWORDS: delta invariant; Tjurina number
1550"
1551{
1552   option(redSB);
1553   def R=basering;
1554   execute("ring S=("+charstr(R)+"),(x,y),lp;");
1555   map phi=R,x,y;
1556   ideal singL=phi(singL);
1557   singL=simplify(std(singL),1);
1558   attrib(singL,"isSB",1);
1559   int d=vdim(singL);
1560   poly f=phi(f);
1561   int i;
1562   int w = printlevel-voice+2;  // w=printlevel (default: w=0)
1563   if(d==1)
1564   {
1565      map alpha=S,var(1)-singL[2][2],var(2)-singL[1][2];
1566      f=alpha(f);
1567      execute("ring C=("+charstr(S)+"),("+varstr(S)+"),ds;");
1568      poly f=imap(S,f);
1569      ideal singL=imap(S,singL);
1570      if((w>=1)&&(ord(f)>=2))
1571      {
1572        "local analysis of the singularities";"";
1573        basering;
1574        singL;
1575        f;
1576        pause();
1577      }
1578   }
1579   else
1580   {
1581      poly p;
1582      poly c;
1583      map psi;
1584      number co;
1585
1586      while((deg(lead(singL[1]))>1)&&(deg(lead(singL[2]))>1))
1587      {
1588         psi=S,x,y+random(-100,100)*x;
1589         singL=psi(singL);
1590         singL=std(singL);
1591          f=psi(f);
1592      }
1593
1594      if(deg(lead(singL[2]))==1)
1595      {
1596         p=singL[1];
1597         c=singL[2]-lead(singL[2]);
1598         co=leadcoef(singL[2]);
1599      }
1600      if(deg(lead(singL[1]))==1)
1601      {
1602         psi=S,y,x;
1603         f=psi(f);
1604         singL=psi(singL);
1605         p=singL[2];
1606         c=singL[1]-lead(singL[1]);
1607         co=leadcoef(singL[1]);
1608      }
1609
1610      execute("ring B=("+charstr(S)+"),a,dp;");
1611      map beta=S,a,a;
1612      poly p=beta(p);
1613
1614      execute("ring C=("+charstr(S)+",a),("+varstr(S)+"),ds;");
1615      number p=number(imap(B,p));
1616      minpoly=p;
1617
1618      map iota=S,a,a;
1619      number c=number(iota(c));
1620      number co=iota(co);
1621
1622      map alpha=S,x-c/co,y+a;
1623      poly f=alpha(f);
1624      f=cleardenom(f);
1625      if((w>=1)&&(ord(f)>=2))
1626      {
1627        "local analysis of the singularities";"";
1628        basering;
1629        alpha;
1630        f;
1631        pause();
1632        "";
1633      }
1634   }
1635   int intMult = deg(lead(f));
1636   poly fdummy = f;
1637   poly gdummy = lead(f);
1638   int ivr = 1;
1639   while(ivr)
1640      {
1641        fdummy = fdummy - lead(fdummy);
1642        if((fdummy ==0) || (deg(lead(fdummy))>intMult)){break;}
1643        gdummy = gdummy + lead(fdummy);
1644      }
1645   ideal SQRpart = sqrfree(gdummy);
1646   int IntForRet;
1647   if(deg(SQRpart[1])==intMult)
1648     {
1649        IntForRet = intMult;
1650     }
1651   option(noredSB);
1652   ideal fstd=std(ideal(f)+jacob(f));
1653   poly hc=highcorner(fstd);
1654   int tau=vdim(fstd);
1655   int o=ord(f);
1656   int delt,nb;
1657
1658   if(tau==0)                 //smooth case
1659   {
1660      setring R;
1661      return(list(0,0,1,0));
1662   }
1663   if((char(basering)>=181)||(char(basering)==0))
1664   {
1665      if(o==2)                //A_k-singularity
1666      {
1667        if(w>=1){"A_k-singularity";"";}
1668         setring R;
1669         delt=(tau+1)/2;
1670         return(list(d*delt,d*tau,d*(2*delt-tau+1),IntForRet));
1671      }
1672      if((lead(f)==var(1)*var(2)^2)||(lead(f)==var(1)^2*var(2)))
1673      {
1674        if(w>=1){"D_k- singularity";"";}
1675
1676         setring R;
1677         delt=(tau+2)/2;
1678         return(list(d*delt,d*tau,d*(2*delt-tau+1),IntForRet));
1679      }
1680
1681      int mu=vdim(std(jacob(f)));
1682      poly g=f+var(1)^mu+var(2)^mu;  //to obtain a convenient Newton-polygon
1683
1684      list NP=newtonpoly(g);
1685      if(w>=1){"Newton-Polygon:";NP;"";}
1686      int s=size(NP);
1687
1688      if(is_NND(f,mu,NP))
1689      { // the Newton-polygon is non-degenerate
1690        // compute nb, the number of branches
1691        for(i=1;i<=s-1;i++)
1692        {
1693          nb=nb+gcd(NP[i][2]-NP[i+1][2],NP[i][1]-NP[i+1][1]);
1694        }
1695        if(w>=1){"Newton-Polygon is non-degenerated";"";}
1696        return(list(d*(mu+nb-1)/2,d*tau,d*nb,IntForRet));
1697      }
1698
1699      if(w>=1){"Newton-Polygon is degenerated";"";}
1700
1701      // the following can certainly be made more efficient when replacing
1702      // 'hnexpansion' (used only for computing number of branches) by
1703      // successive blowing-up + test if Newton polygon degenerate:
1704      if(s>2)    //  splitting of f
1705      {
1706         if(w>=1){"Newton polygon can be used for splitting";"";}
1707         intvec v=NP[1][2]-NP[2][2],NP[2][1];
1708         int de=w_deg(g,v);
1709         int st=w_deg(hc,v)+v[1]+v[2];
1710         poly f1=var(2)^NP[2][2];
1711         poly f2=jet(g,de,v)/var(2)^NP[2][2];
1712         poly h=g-f1*f2;
1713         de=w_deg(h,v);
1714         poly k;
1715         ideal wi=var(2)^NP[2][2],f2;
1716         matrix li;
1717         while(de<st)
1718         {
1719           k=jet(h,de,v);
1720           li=lift(wi,k);
1721           f1=f1+li[2,1];
1722           f2=f2+li[1,1];
1723           h=g-f1*f2;
1724           de=w_deg(h,v);
1725         }
1726         nb=deltaLocMod(f1,maxideal(1))[3]+deltaLocMod(f2,maxideal(1))[3];
1727         setring R;
1728         return(list(d*(mu+nb-1)/2,d*tau,d*nb,IntForRet));
1729      }
1730
1731      f=jet(f,deg(hc)+2);
1732      if(w>=1){"now we have to use Hamburger-Noether (Puiseux) expansion";}
1733      ideal fac=factorize(f,1);
1734      if(size(fac)>1)
1735      {
1736         nb=0;
1737         for(i=1;i<=size(fac);i++)
1738         {
1739            nb=nb+deltaLocMod(fac[i],maxideal(1))[3];
1740         }
1741         setring R;
1742         return(list(d*(mu+nb-1)/2,d*tau,d*nb,IntForRet));
1743      }
1744      list HNEXP=hnexpansion(f);
1745      if (typeof(HNEXP[1])=="ring") {
1746        def altring = basering;
1747        def HNEring = HNEXP[1]; setring HNEring;
1748        nb=size(hne);
1749        setring R;
1750        kill HNEring;
1751      }
1752      else
1753      {
1754        nb=size(HNEXP);
1755      }
1756      return(list(d*(mu+nb-1)/2,d*tau,d*nb,IntForRet));
1757   }
1758   else             //the case of small characteristic
1759   {
1760      f=jet(f,deg(hc)+2);
1761      if(w>=1){"now we have to use Hamburger-Noether (Puiseux) expansion";}
1762      delt=delta(f);
1763      return(list(d*delt,d*tau,d,IntForRet));
1764   }
1765}
1766
1767///////////////////////////////////////////////////////////////////////////////
1768proc paraConic(poly q)
1769"USAGE:  paraConic(q); q poly
1770ASSUME:  The basering must be a polynomial ring in three variables with
1771         coefficients in Q. @*
1772         The polynomial q must be homogeneous of degree 2 and absolutely
1773         irreducible. @*
1774NOTE:    The procedure might fail or give a wrong output if the assumptions
1775         do not hold.
1776
1777RETURN:  ring with an ideal PARACONIC. The ring should be considered as the
1778         homogeneous coordinate ring of PP^1, the ideal defines a rational
1779         parametrization PP^1 --> C2 = {q=0}.
1780
1781THEORY:  We compute a point on C2 via @ref{rationalPointConic}. The pencil of
1782         lines through this point projects C2 birationally to PP^1. Inverting
1783         the projection gives the result.
1784KEYWORDS: conic, parametrization, rational point.
1785SEE ALSO: rationalPointConic.
1786EXAMPLE: example paraConic; shows an example
1787"
1788{
1789  def Roriginal = basering;
1790  def RP2 = projConic(q);  // ring with ideal Ipoint
1791                           // possibly defined over algebraic number field
1792  setring RP2;
1793  def rp1 = invertBirMap(Ipoint, ideal(fetch(Roriginal,q)));
1794  setring rp1;
1795  list rl = ringlist(rp1);
1796  rl[2] = list("s","t");
1797  def RP1 = ring(rl);
1798  setring RP1;
1799  ideal PARACONIC = fetch(rp1,psi);
1800  export(PARACONIC);
1801  "// 'paraConic' created a ring together with an ideal RNC.";
1802  "// Supposing you typed, say,  def RP1 = paraConic(q);";
1803  "// you may access the ideal by typing";
1804  "//      setring RP1; PARACONIC;";
1805  return(RP1);
1806}
1807example
1808{ "EXAMPLE:"; echo=2;
1809  ring R = 0,(x,y,z),dp;
1810  poly f = y^8-x^3*(z+x)^5;
1811  ideal adj = adjointIdeal(f);
1812  def Rn = invertBirMap(adj,ideal(f));
1813  setring(Rn);
1814  J;
1815  def Rc = rncItProjEven(J);
1816  PHI;
1817  setring Rc;
1818  CONIC;
1819  def RPc = paraConic(CONIC);
1820  setring RPc;
1821  PARACONIC;
1822}
1823
1824///////////////////////////////////////////////////////////////////////////////
1825static proc projConic(poly q)
1826"USAGE:  projConic(q); q poly
1827ASSUME:  The basering must be a polynomial ring in three variables with
1828         coefficients in Q. @*
1829         The polynomial q must be homogeneous of degree 2 and absolutely
1830         irreducible. @*
1831NOTE:    The procedure might fail or give a wrong output if the assumptions
1832         do not hold.
1833RETURN:  ring with an ideal Ipoint defining a pencil of lines through a point
1834         on the conic C2 = {q=0}. This point has either coefficients in Q or
1835         in a quadratic extension field of Q.
1836THEORY:  We compute the point on C2 via @ref rationalPointConic.
1837KEYWORDS: conic, parametrization, rational point.
1838SEE ALSO: rationalPointConic.
1839"
1840{
1841  def Roriginal = basering;
1842  list rl = ringlist(Roriginal);
1843  rl[3] = list(list("dp",1:3),list("C",0));
1844  def RP20 = ring(rl);
1845  setring RP20;
1846  poly q = imap(Roriginal,q);
1847  def RP21 = rationalPointConic(q);  //  ring with ideal point representing
1848                                     //  point on conic
1849                                     //  possibly defined over algebraic number
1850                                     //  field
1851  setring RP21;
1852  list rl1 = ringlist(RP21);
1853  rl1[2] = list("u","v","w");
1854  rl1[3] = list(list("dp",1:3),list("C",0));
1855  def RP2 = ring(rl1);
1856  setring RP2;
1857  ideal point = fetch(RP21,point);
1858  matrix bP = syz(point);
1859  ideal Ipoint = matrix(maxideal(1))*bP;  // defines pencil of lines through
1860                                          // point
1861  export(Ipoint);
1862  return(RP2);
1863}
1864
1865
1866///////////////////////////////////////////////////////////////////////////////
1867static proc isIrreducible(poly f)
1868"USAGE:  isIrreducible(f); f poly
1869RETURN:  1 iff the given polynomial f is irreducible; 0 otherwise.
1870THEORY:  This test is performed by computing the absolute factorization of f.
1871KEYWORDS: irreducible.
1872"
1873{
1874  def r = basering;
1875  def s = absFactorize(f);
1876  setring s;
1877  list L = absolute_factors;
1878  int result = 0;
1879  if (L[4] == 1){result = 1;}
1880  setring r;
1881  kill s;
1882  return (result);
1883}
1884
1885
1886///////////////////////////////////////////////////////////////////////////////
1887static proc isQuadratic(poly p)
1888"USAGE:  isQuadratic(p); p poly
1889RETURN:  checks whether p is a homogeneous, quadratic polynomial in the
1890         first three ring variables, x, y, z say;
1891         If so, the method extracs the coefficients a, b, ..., f such that
1892         p = a*x2 + b*xy + c * y2 + d * xz + e * yz + f * z2
1893         and returns them as a list of seven entries, [1, a, b, c, d, e, f];
1894         otherwise, a list with the single entry [0] is returned
1895"
1896{
1897  bigint a = bigint(leadcoef(subst(p, var(2), 0, var(3), 0)));
1898  bigint c = bigint(leadcoef(subst(p, var(1), 0, var(3), 0)));
1899  bigint f = bigint(leadcoef(subst(p, var(1), 0, var(2), 0)));
1900  poly h = p - a * var(1)^2 - c * var(2)^2 - f * var(3)^2;
1901  bigint b = bigint(leadcoef(subst(h, var(3), 0)));
1902  bigint d = bigint(leadcoef(subst(h, var(2), 0)));
1903  bigint e = bigint(leadcoef(subst(h, var(1), 0)));
1904  list L = 0;
1905  if (h - b * var(1) * var(2) - d * var(1) * var(3)
1906        - e * var(2) * var(3) != 0) { return (L); }
1907  L = 1, a, b, c, d, e, f;
1908  return (L);
1909}
1910
1911
1912///////////////////////////////////////////////////////////////////////////////
1913static proc largestSquare(bigint n)
1914"USAGE:  largestSquare(n); n bigint
1915ASSUME:  n <> 0
1916RETURN:  returns the largest positive number m (as bigint) such that m^2
1917         divides n.
1918THEORY:  This computation is done by prime factorization of n.
1919KEYWORDS: prime factorization.
1920"
1921{
1922  if (n == 0) { "ERROR: largestSquare(0) had been invoked"; }
1923
1924  bigint nn = n; if (nn < 0) { nn = -n; }
1925  list L = primefactors(nn);
1926  if (L[3] != 1)
1927  { "WARNING: command 'primefactors(.)' did not find all prime factors"; }
1928  int i; bigint m = bigint(1); int e; int j;
1929  for (i = 1; i <= size(L[1]); i++)
1930  {
1931    e = L[2][i] / 2;
1932    for (j = 1; j <= e; j++) { m = m * bigint(L[1][i]); }
1933  }
1934  return (m);
1935}
1936
1937
1938///////////////////////////////////////////////////////////////////////////////
1939static proc jIndex(bigint a, bigint b, bigint c)
1940"USAGE:  jIndex(a, b, c); a, b, c bigint's
1941RETURN:  returns the middle of the three numbers |ab|, |bc|, and |ca|.
1942"
1943{
1944  bigint n1 = a*b; if (n1 < 0) { n1 = -n1; }
1945  bigint n2 = b*c; if (n2 < 0) { n2 = -n2; }
1946  bigint n3 = c*a; if (n3 < 0) { n3 = -n3; }
1947  if ((n1 <= n2) && (n2 <= n3)) { return (n2); }
1948  if ((n1 >= n2) && (n2 >= n3)) { return (n2); }
1949  if ((n2 <= n1) && (n1 <= n3)) { return (n1); }
1950  if ((n2 >= n1) && (n1 >= n3)) { return (n1); }
1951  return (n3);
1952}
1953
1954
1955///////////////////////////////////////////////////////////////////////////////
1956static proc polyModP(poly q, bigint p)
1957"USAGE:  polyModP(q, p); q poly, p bigint
1958RETURN:  takes each coefficient of q modulo p and returns the resulting poly
1959"
1960{
1961  poly qq = q; poly res = 0;
1962  bigint c;
1963  while (qq != 0)
1964  {
1965    c = bigint(leadcoef(qq) mod p);
1966    res = res + c * leadmonom(qq);
1967    qq = qq - lead(qq);
1968  }
1969  return (res);
1970}
1971
1972
1973///////////////////////////////////////////////////////////////////////////////
1974static proc rootModP(bigint r, bigint p)
1975"USAGE:  rootModP(r, p); r, p bigint's
1976ASSUME:  0 <= r < p, and p prime;
1977         Furthermore it is assumes that there is some x in {0, 1, ..., p-1}
1978         such that x^2 = r mod p;
1979RETURN:  an x in {0, 1, ..., p-1} such that x^2 = r mod p;
1980THEORY:  For p larger than 32003, this computation is done using Cantor-
1981         Zassenhaus' algorithm. Otherwise a brute force approach is used.
1982KEYWORDS: Cantor-Zassenhaus algorithm.
1983"
1984{
1985  if (r == 0) { return (0); }
1986  if (r == 1) { return (1); }
1987  if (p <= 32003)
1988  {
1989    /* For small p, we use a brute force approach: */
1990    int i;
1991    for (i = 2; i < p; i++)
1992    {
1993      if (((i*i) mod p) == r) { return (i); }
1994    }
1995    /* should never be reached: */
1996    return (-1);
1997  }
1998
1999  /* For p > 32003, we use Cantor-Zassenhaus' algorithm: */
2000  def br = basering;
2001  ring rTemp = 0, x, dp;
2002  bigint b; bigint exponent; poly factor;
2003  poly h = x^2 - r;
2004  ideal redI = h; redI = std(redI);
2005  poly q = x^2; bigint root = 0;
2006  while (root == 0)
2007  {
2008    b = bigint(random(1, 2^30));
2009    exponent = bigint((p - 1) / 2);
2010    /* We need to compute q^exponent mod (x^2 - a) and mod p: */
2011    factor = x + b; q = 1;
2012    while (exponent > 0)
2013    {
2014      if ((exponent mod 2) == 1)
2015      {
2016        q = q * factor;
2017        q = reduce(q, redI);
2018        q = polyModP(q, p);
2019      }
2020      exponent = bigint(exponent / 2);
2021      factor = factor * factor;
2022      factor = reduce(factor, redI);
2023      factor = polyModP(factor, p);
2024    }
2025    if (deg(q) == 1)
2026    {
2027      q = q - 1;
2028      b = inverseModP(bigint(leadcoef(q)), p);
2029      q = q - lead(q);
2030      root = (bigint(q) * b) mod p;
2031      if (((root * root - r) mod p) != 0) { root = 0; }
2032    }
2033  }
2034  setring br; kill rTemp;
2035  return (root);
2036}
2037
2038
2039///////////////////////////////////////////////////////////////////////////////
2040static proc inverseModP(bigint r, bigint p)
2041"USAGE:  inverseModP(r, p); r, p bigint's
2042ASSUME:  0 <= r < p, and r and p coprime;
2043RETURN:  returns the inverse of r in Z/p represented by an element in
2044         {1, 2, ..., p-1}
2045THEORY:  This uses Euclid's extended gcd algorithm.
2046"
2047{
2048  list L = extendedEuclid(r, p);
2049  if (L[1] != 1) { "ERROR: GCD of", r, "and", p, "should be 1."; }
2050  L[2] = L[2] mod p;
2051  return (L[2]);
2052}
2053
2054
2055///////////////////////////////////////////////////////////////////////////////
2056static proc squareRoot(bigint r, bigint m, int justCheck)
2057"USAGE:  squareRoot(r, m, j); r, m bigint's, j int
2058RETURN:  checks whether r is a square modulo m, i.e., checks whether there is
2059         some x such that x^2 = r mod m;
2060         If justCheck is 1, then the method will terminate after the check
2061         and return 1 if r is a square and -1 otherwise.
2062         If justCheck is 0 and r is a square, then the method continues and
2063         computes a solution x in {0, 1, m-1} with x^2 = r mod m, which will
2064         then be returned
2065THEORY:  This algorithm checks solvability by computing the Legendre symbols
2066         modulo all primes in m. Then, individual roots will be computed and
2067         lifted to the desired square root modulo m using Chinese
2068         remaindering.
2069"
2070{
2071  if (m == 0) { "ERROR: squareRoot had been invoked with m = 0"; }
2072
2073  list L = primefactors(m);
2074  if ((L[3] != 1) && (L[3] != -1))
2075  { "WARNING: command 'primefactors(.)' did not find all prime factors"; }
2076  int i;
2077  for (i = 1; i <= size(L[2]); i++)
2078  {
2079    if (legendreSymbol(r, L[1][i]) == -1) { return (-1); }
2080  }
2081  /* now we know that there is some x in {0, 1, m-1} with
2082     x^2 = r mod m */
2083  if (justCheck == 1) { return (1); }
2084  else
2085  {
2086    // now we need to compute x; this works in two stages:
2087    // 1) write m = p1^e1 * ... * pk^ek (we already have that),
2088    // 2) for each i in {1, 2, ..., k}
2089    //    2.1) compute a yi such that yi^2 = r mod pi,
2090    //    2.2) lift yi to an xi such that xi^2 = r mod (pi^ei),
2091    // 3) lift (x1, x2, ..., xk) in Z/p1^e1 * ... * Z/pk^ek
2092    //    to x in Z/m via Chinese remainder theorem
2093
2094    list roots;
2095    // 2.1):
2096    for (i = 1; i <= size(L[1]); i++)
2097    {
2098      roots = insert(roots, rootModP(r mod L[1][i], L[1][i]), size(roots));
2099    }
2100
2101    // 2.2):
2102    bigint c; bigint l; bigint temp; bigint pPower; int e;
2103    for (i = 1; i <= size(roots); i++)
2104    {
2105      pPower = bigint(L[1][i]);
2106      for (e = 2; e <= L[2][i]; e++)
2107      {
2108        c = bigint(roots[i]); l = pPower;
2109        temp = r - c * c; l = bigint(2) * c * l; c = temp;
2110        c = c div pPower; l = l div pPower;
2111        c = c mod L[1][i]; l = l mod L[1][i];
2112        c = (c * bigint(inverseModP(l, L[1][i]))) mod L[1][i];
2113        c = bigint(roots[i]) + c * pPower;
2114        pPower = pPower * L[1][i]; roots[i] = c;
2115      }
2116    }
2117
2118    // 2.3):
2119    list mm; bigint z; int j;
2120    for (i = 1; i <= size(L[1]); i++)
2121    {
2122      z = bigint(L[1][i]);
2123      for (j = 2; j <= L[2][i]; j++)
2124      {
2125        z = z * bigint(L[1][i]);
2126      }
2127      mm = insert(mm, z, size(mm));
2128    }
2129    return (chineseRemainder(roots, mm) mod m);
2130  }
2131}
2132
2133
2134///////////////////////////////////////////////////////////////////////////////
2135static proc chineseRemainder(list rr, list mm)
2136"USAGE:  chineseRemainder(rr, mm); rr, mm lists of bigint's
2137ASSUME:  lists rr and mm must have same sizes;
2138         Furthermore the entries of mm must be mutually coprime.
2139RETURN:  an x which fulfills the simultaneous remainder conditions
2140         x = rr[i] mod mm[i], 1 <= i <= size(rr)
2141KEYWORDS: Chinese remainder.
2142"
2143{
2144  bigint x = bigint(0); int i; bigint N; list l;
2145  bigint M = bigint(mm[1]);
2146  for (i = 2; i <= size(mm); i++) { M = M * bigint(mm[i]); }
2147  for (i = 1; i <= size(mm); i++)
2148  {
2149    N = M div mm[i];
2150    l = extendedEuclid(mm[i], N);
2151    x = x + rr[i]*l[3]*N;
2152  }
2153  return (x);
2154}
2155
2156
2157///////////////////////////////////////////////////////////////////////////////
2158static proc rationalPointSpecial(bigint b1, bigint c1)
2159"USAGE:  rationalPointSpecial(b1, c1); b1, c1 bigint's
2160ASSUME:  b1 <> 0 and c1 <> 0;
2161RETURN:  with poly p = var(1)^2 + b1 * var(2)^2 + c1 * var(3)^2, the method
2162         returns a list L with either one entry or four entries:
2163         case 'three entries':
2164           L[1] = 0 signaling that there is no rational point on V(p),
2165           L[2] the largest number b such that b^2 divides b1
2166                (for subsequent use by the caller of this method),
2167           L[3] the largest number c such that c^2 divides c1
2168                (for subsequent use by the caller of this method);
2169         case 'four entries':
2170           L[1] = 1 signaling that there is a rational point on V(p),
2171           L[2], L[3], L[4] rational numbers such that the tuple
2172                (L[2], L[3], L[4]) is on V(p)
2173"
2174{
2175  if (b1 == 0) { "ERROR: rationalPointSpecial(0, c1) had been invoked"; }
2176  if (c1 == 0) { "ERROR: rationalPointSpecial(b1, 0) had been invoked"; }
2177
2178  bigint b_s = largestSquare(b1); bigint b_r = b1/b_s/b_s;
2179  bigint c_s = largestSquare(c1); bigint c_r = c1/c_s/c_s;
2180  bigint g = gcd(b_r, c_r);
2181  def S=basering;
2182  ideal mi = maxideal(1);
2183  map mm = basering, mi; map mTemp;
2184  mm[1] = var(1); mm[2] = var(2)/b_s/g; mm[3] = var(3)/c_s/g;
2185  bigint a = g;     bigint aa = a; if (aa <= 0) { aa = -aa; }
2186  bigint b = b_r/g; bigint bb = b; if (bb <= 0) { bb = -bb; }
2187  bigint c = c_r/g; bigint cc = c; if (cc <= 0) { cc = -cc; }
2188  bigint R1 = squareRoot(-a*b, cc, 1);
2189  if (R1 == -1) { list L = 0, b_s, c_s; return (L); }
2190  bigint R2 = squareRoot(-a*c, bb, 1);
2191  if (R2 == -1) { list L = 0, b_s, c_s; return (L); }
2192  bigint R3 = squareRoot(-b*c, aa, 1);
2193  if (R3 == -1) { list L = 0, b_s, c_s; return (L); }
2194  bigint t; bigint r1; bigint Q; bigint A; bigint B; bigint C;
2195  bigint alpha; bigint beta; bigint gamma;
2196  while (jIndex(a, b, c) > 1)
2197  {
2198    mTemp = basering, mi;
2199    if (aa > cc)
2200    {
2201      t = a; a = c; c = t;
2202      t = aa; aa = cc; cc = t;
2203      mTemp = basering, mi;
2204      mTemp[1] = var(3); mTemp[3] = var(1); mm = mTemp(mm);
2205    }
2206    if (bb > cc)
2207    {
2208      t = b; b = c; c = t;
2209      t = bb; bb = cc; cc = t;
2210      mTemp = basering, mi;
2211      mTemp[2] = var(3); mTemp[3] = var(2); mm = mTemp(mm);
2212    }
2213    if (bb < aa)
2214    {
2215      t = b; b = a; a = t;
2216      t = bb; bb = aa; aa = t;
2217      mTemp = basering, mi;
2218      mTemp[1] = var(2); mTemp[2] = var(1); mm = mTemp(mm);
2219    }
2220    /* now, we have established |a| <= |b| <= |c|; and permuted
2221       the map mm, accordingly */
2222    cc = c; if (cc <= 0) { cc = -cc; }
2223    R1 = squareRoot(-a*b, cc, 0);
2224    r1 = (R1 * inverseModP(a, cc)) mod cc;
2225    if (r1*bigint(2) > cc) { r1 = r1 - cc; }
2226    Q = (a*r1*r1 + b)/c;
2227    if (Q == 0)
2228    {
2229      list L = 1, subst(mm[1], var(1), 1, var(2), 1, var(3), 0),
2230                  subst(mm[2], var(1), 1, var(2), 1, var(3), 0),
2231                  subst(mm[3], var(1), 1, var(2), 1, var(3), 0);
2232      return (L);
2233    }
2234    A = gcd(gcd(a*r1*r1, b), c*Q);
2235    alpha = r1/A; beta = b/A;
2236    B = a*beta;
2237    gamma = largestSquare(Q/A);
2238    C = Q/A/gamma/gamma;
2239    mTemp = basering, mi;
2240    mTemp[1] = A*alpha*var(1) - beta*var(2);
2241    mTemp[2] = var(1) + a*alpha*var(2);
2242    mTemp[3] = C*gamma*var(3);
2243    mm = mTemp(mm);
2244    a = A; b = B; c = C;
2245    aa = a; if (aa <= 0) { aa = -aa; }
2246    bb = b; if (bb <= 0) { bb = -bb; }
2247    cc = c; if (cc <= 0) { cc = -cc; }
2248  }
2249  if (a*b < 0)
2250  {
2251    list L = 1, subst(mm[1], var(1), 1, var(2), 1, var(3), 0),
2252                subst(mm[2], var(1), 1, var(2), 1, var(3), 0),
2253                subst(mm[3], var(1), 1, var(2), 1, var(3), 0);
2254    return (L);
2255  }
2256  if (a*c < 0)
2257  {
2258    list L = 1, subst(mm[1], var(1), 1, var(2), 0, var(3), 1),
2259                subst(mm[2], var(1), 1, var(2), 0, var(3), 1),
2260                subst(mm[3], var(1), 1, var(2), 0, var(3), 1);
2261    return (L);
2262  }
2263  if (b*c < 0)
2264  {
2265    list L = 1, subst(mm[1], var(1), 0, var(2), 1, var(3), 1),
2266                subst(mm[2], var(1), 0, var(2), 1, var(3), 1),
2267                subst(mm[3], var(1), 0, var(2), 1, var(3), 1);
2268    return (L);
2269  }
2270  list L = 0, b_s, c_s; return (L);
2271}
2272
2273
2274///////////////////////////////////////////////////////////////////////////////
2275static proc extendedEuclid(bigint a, bigint b)
2276"USAGE:  extendedEuclid(a, b); a, b bigint's
2277ASSUME:  a <> 0 or b <> 0;
2278RETURN:  returns a list with three entries:
2279         _[1]: gcd(a,b) > 0,
2280         _[2], _[3]: s, t, such that s*a + t*b = gcd(a,b)
2281KEYWORDS: extended Euclidean algorithm.
2282"
2283{
2284  list l = 0; bigint temp;
2285  if (a == 0)         { l = b, 0, 1; if (b < 0) { l = -b, 0, -1; } }
2286  if (b == 0)         { l = a, 1, 0; if (a < 0) { l = -a, -1, 0; } }
2287  if ((a mod b) == 0) { l = b, 0, 1; if (b < 0) { l = -b, 0, -1; } }
2288  if ((b mod a) == 0) { l = a, 1, 0; if (a < 0) { l = -a, -1, 0; } }
2289  if (size(l) > 1) { return (l); }
2290
2291  temp = a mod b;
2292  l = extendedEuclid(b, temp);
2293  temp = (a - temp) / b;
2294  temp = bigint(l[2]) - temp * bigint(l[3]);
2295  l = l[1], l[3], temp;
2296  return (l);
2297}
2298
2299static proc legendreSymbol(bigint r, bigint p)
2300"assumes p prime;
2301returns the Legendre symbol (r/p), that is
2302 1 if r appears as residue modulo p of a square,
2303-1 if not,
2304 0 if r is a multiple of p
2305"
2306{
2307  bigint rr = r mod p;
2308  if (rr == 0) { return (0) }
2309  if (rr == 1) { return (1) }
2310  /* now, p must be at least 3 */
2311  bigint e = (p - 1) / bigint(2);
2312  bigint result = 1;
2313  bigint power = rr;
2314  while (e > 0)
2315  {
2316    if ((e mod 2) == 1) { result = (result * power) mod p; }
2317    e = e / bigint(2);
2318    power = (power * power) mod p;
2319  }
2320  if (result > 1) { result = result - p; /* should be -1 */ }
2321  return (result);
2322}
2323
2324
2325///////////////////////////////////////////////////////////////////////////////
2326static proc buildExtension(bigint b, bigint c, bigint bs, bigint cs)
2327"USAGE:  buildExtension(b, c, bs, cs); b, c, bs, cs bigint's
2328ASSUME:  assumes that bs is the largest positive number such that bs^2
2329         divides b; analogously for cs regarding c;
2330         Assumes furthermore that there is no rational point on the conic
2331         X^2 + b*Y^2 + c*Z^2 = 0.
2332         Assumes that the ground field of the basering is Q.
2333RETURN:  builds an appropriate quadratic field extension Q(a) in which a
2334         point exists that lies on the given conic. This point is stored in
2335         a newly defined and exported (1x3) matrix named 'point'.
2336         The method returns the resulting polynomial ring over Q(a).
2337"
2338{
2339  bigint br = b/bs/bs;
2340  bigint cr = c/cs/cs;
2341  /* X^2 + br*bs^2*Y^2 + cr*cs^2*Z^2 = 0 */
2342  def bRing = basering;
2343  list L = ringlist(bRing);
2344
2345  if (b != 0)
2346  {
2347    L[1] = list(0, list("a"), list(list("lp", 1)), ideal(0));
2348    def RTemp = ring(L);
2349    setring RTemp; list L = ringlist(RTemp);
2350    L[1][4] = ideal(a^2 + br);
2351    def R = ring(L);
2352    setring R; kill RTemp;
2353    matrix point[1][3];
2354    point[1, 1] = a * bs; point[1, 2] = 1; point[1, 3] = 0;
2355    export point;
2356    setring bRing;
2357    return (R);
2358  }
2359  if (c != 0)
2360  {
2361    L[1] = list(0, list("a"), list(list("lp", 1)), ideal(0));
2362    def RTemp = ring(L);
2363    setring RTemp; list L = ringlist(RTemp);
2364    L[1][4] = ideal(a^2 + cr);
2365    def R = ring(L);
2366    setring R; kill RTemp;
2367    matrix point[1][3];
2368    point[1, 1] = a * cs; point[1, 2] = 0; point[1, 3] = 1;
2369    export point;
2370    setring bRing;
2371    return (R);
2372  }
2373
2374  "ERROR: unexpectedly encountered conic X^2 + 0*Y^2 + 0*Z^2 = 0";
2375  return (bRing);
2376}
2377
2378
2379///////////////////////////////////////////////////////////////////////////////
2380static proc testRationalPointConic(poly pp)
2381"USAGE:  testRationalPointConic(pp); pp poly
2382RETURN:  returns 0 in case of unexpected input (e.g. non-quadratic,
2383         reducible); 1 otherwise
2384NOTE:    This method calles rationalPointConic, measures time consumption
2385         and checks whether the computed point lies indeed on the conic pp.
2386         The results will be printed to standard output.
2387"
2388{
2389  "testing rationalPointConic(poly) for poly p:";
2390  "p =", pp;
2391  if (isQuadratic(pp)[1] == 1) { "p is quadratic."; }
2392  else                         { "p is not quadratic.";   return (0); }
2393  if (isIrreducible(pp) == 1)  { "p is irreducible."; }
2394  else                         { "p is not irreducible."; return (0); }
2395  int tps = system("--ticks-per-sec");
2396  system("--ticks-per-sec", 1000);
2397  def rOrig = basering;
2398  int t = rtimer;
2399  def rNew = rationalPointConic(pp);
2400  t = rtimer - t;
2401  "time for finding a point on the conic [msec] =", t;
2402  system("--ticks-per-sec", tps);
2403  setring rNew;
2404  poly ff = fetch(rOrig, pp);
2405  if (minpoly == 0)
2406  { "there is a rational point on the conic p";
2407    "x =", point[1,1], "  y =", point[1,2], "  z =", point[1,3];
2408    "check (should be zero):", subst(ff, var(1), point[1,1],
2409                                         var(2), point[1,2],
2410                                         var(3), point[1,3]);
2411  }
2412  else
2413  {
2414    "there is no rational point on the conic p";
2415    "but there is a point on the conic in the field extension Q(a),";
2416    "with minpoly =", minpoly;
2417    "x =", point[1,1], "  y =", point[1,2], "  z =", point[1,3];
2418    "check (should be zero):", subst(ff, var(1), point[1,1],
2419                                         var(2), point[1,2],
2420                                         var(3), point[1,3]);
2421  }
2422  setring rOrig;
2423}
2424
2425example
2426{ "EXAMPLE:"; echo=2;
2427  ring r = 0, (x,y,z, u, v, w), dp;
2428  poly p = x^2 + 2*y^2 + 5*z^2 - 4*x*y + 3*x*z + 17*y*z;
2429  testRationalPointConic(p);
2430}
2431
2432///////////////////////////////////////////////////////////////////////////////
2433proc rationalPointConic(poly p)
2434"USAGE:  rationalPointConic(p); p poly
2435ASSUME:  assumes that p is an irreducible quadratic polynomial in the first
2436         three ring variables;
2437         ground field is expected to be Q.
2438RETURN:  The method finds a point on the given conic. There are two
2439         possibilities:
2440         1) There is a rational point on the curve.
2441         2) There is no rational point on the curve.
2442         In the second case, the method creates a modification of the current
2443         basering which is a polynomial ring over some quadratic field
2444         extension Q(a) of Q. Apart from the replacement of Q by Q(a), the
2445         new polynomial ring, R say, is the same as the original basering.
2446         (In the first case, R is identical with the basering.)
2447         In both cases, the method will then define a (1x3) matrix named
2448         'point' which lives in R and which contains the coordinates of the
2449         desired point on q.
2450         Finally, the method returns the ring R (which will in the 1st case
2451         be the original base ring).
2452EXAMPLE: example rationalPointConic; shows an example
2453"
2454{
2455  list L = isQuadratic(p);
2456  bigint a = bigint(L[2]); bigint b = bigint(L[3]); bigint c = bigint(L[4]);
2457  bigint d = bigint(L[5]); bigint e = bigint(L[6]); bigint f = bigint(L[7]);
2458  bigint x; bigint y; bigint z; bigint nn;
2459  def R = basering;
2460
2461  if (b^2 == 4*a*c)
2462  {
2463    if (c == 0)
2464    {
2465      x = -2*d*e; y = d^2-4*a*f; z = e*4*a;
2466      nn = gcd(gcd(absValue(x), absValue(y)), absValue(z));
2467      matrix point[1][3];
2468      point[1, 1] = x/nn; point[1, 2] = y/nn; point[1, 3] = z/nn;
2469      export point; return (R);
2470    }
2471    else
2472    {
2473      bigint fs = 4*c*f - e^2;
2474      bigint ds = 4*c*d - 2*b*e;
2475      x = -fs*2*c; y = b*fs-e*ds; z = ds*2*c;
2476      nn = gcd(gcd(absValue(x), absValue(y)), absValue(z));
2477      matrix point[1][3];
2478      point[1, 1] = x/nn; point[1, 2] = y/nn; point[1, 3] = z/nn;
2479      export point; return (R);
2480    }
2481  }
2482
2483  if (d^2 == 4*a*f)
2484  {
2485    if (f == 0)
2486    {
2487      x = -b*e*2; y = e*4*a; z = b^2-4*a*c;
2488      nn = gcd(gcd(absValue(x), absValue(y)), absValue(z));
2489      matrix point[1][3];
2490      point[1, 1] = x/nn; point[1, 2] = y/nn; point[1, 3] = z/nn;
2491      export point; return (R);
2492    }
2493    else
2494    {
2495      bigint c_s = 4*c*f - e^2;
2496      bigint b_s = 4*f*b - 2*d*e;
2497      x = -c_s*2*f; y = b_s*2*f; z = d*c_s-e*b_s;
2498      nn = gcd(gcd(absValue(x), absValue(y)), absValue(z));
2499      matrix point[1][3];
2500      point[1, 1] = x/nn; point[1, 2] = y/nn; point[1, 3] = z/nn;
2501      export point; return (R);
2502    }
2503  }
2504
2505  if (e^2 == 4*c*f)
2506  {
2507    if (c == 0)
2508    {
2509      x = b*4*f; y = d^2-4*a*f; z = -b*d*2;
2510      nn = gcd(gcd(absValue(x), absValue(y)), absValue(z));
2511      matrix point[1][3];
2512      point[1, 1] = x/nn; point[1, 2] = y/nn; point[1, 3] = z/nn;
2513      export point; return (R);
2514    }
2515    else
2516    {
2517      bigint as = 4*c*a - b^2;
2518      bigint ds = 4*c*d - 2*b*e;
2519      x = ds*2*c; y = e*as-b*ds; z = -as*2*c;
2520      nn = gcd(gcd(absValue(x), absValue(y)), absValue(z));
2521      matrix point[1][3];
2522      point[1, 1] = x/nn; point[1, 2] = y/nn; point[1, 3] = z/nn;
2523      export point; return (R);
2524    }
2525  }
2526
2527  ideal mi = maxideal(1);
2528  map mm = R, mi;
2529  bigint B; bigint C; bigint D;
2530
2531  if ((a == 0) && (c == 0))
2532  {
2533    B = -1; C = 4*b*f - 4*d*e;
2534    /* now, b <> 0 since otherwise p would have the factor z,
2535       and hence not be irreducible */
2536    mm[1] = (var(1)+var(2)-2*e*var(3))/(2*b);
2537    mm[2] = (var(1)-var(2)-2*d*var(3))/(2*b);
2538  }
2539  if ((a != 0) && (c == 0))
2540  {
2541    mm[1] = var(2);
2542    mm[2] = var(1);
2543    bigint t = a; a = c; c = t;
2544    t = e; e = d; d = t;
2545  }
2546  if (c != 0)
2547  {
2548    D = 4*a*c-b^2;
2549    mm[2] = (var(2)-e*var(3)-b*var(1))/(2*c);
2550    map mTemp = basering, mi;
2551    mTemp[1] = (var(1)-2*d*c*var(3)+b*e*var(3))/D;
2552    mm = mTemp(mm);
2553    B = D;
2554    C = 16*a*c^2*f-4*a*c*e^2-4*b^2*c*f+4*b*c*d*e-4*c^2*d^2;
2555  }
2556  list K;
2557  if ((B > 0) && (C >= 0)) { K = 0; }
2558  if ((B >= 0) && (C > 0)) { K = 0; }
2559  if (B == 0)
2560  {
2561    /* looking for a point on X^2 = |C| * Z^2 */
2562    bigint root = largestSquare(absValue(C));
2563    if (absValue(C)/root/root == 1) { K = 1, root, 0, 1; }
2564    else                            { K = 0; }
2565  }
2566  if (C == 0)
2567  {
2568    /* looking for a point on X^2 = |B| * Y^2 */
2569    bigint root = largestSquare(absValue(B));
2570    if (absValue(B)/root/root == 1) { K = 1, root, 1, 0; }
2571    else                            { K = 0; }
2572  }
2573  else { K = rationalPointSpecial(B, C); }
2574  if (K[1] == 0)
2575  {
2576    /* no rational point on conic;
2577       we need to move to an appropriate field extension Q(a) */
2578    poly h1 = mm[1]; poly h2 = mm[2]; poly h3 = mm[3];
2579    def extendedR = buildExtension(B, C, K[2], K[3]);
2580    setring extendedR;
2581    poly g1 = fetch(R, h1);
2582    poly g2 = fetch(R, h2);
2583    poly g3 = fetch(R, h3);
2584    matrix temp[1][3];
2585    temp[1, 1] = subst(g1, var(1), point[1, 1], var(2), point[1, 2],
2586                           var(3), point[1, 3]);
2587    temp[1, 2] = subst(g2, var(1), point[1, 1], var(2), point[1, 2],
2588                           var(3), point[1, 3]);
2589    temp[1, 3] = subst(g3, var(1), point[1, 1], var(2), point[1, 2],
2590                           var(3), point[1, 3]);
2591    point[1, 1] = temp[1, 1]; point[1, 2] = temp[1, 2];
2592    point[1, 3] = temp[1, 3];
2593    setring R;
2594    return (extendedR);
2595  }
2596  else
2597  {
2598    string dummyString = string(K); // without this useless line, we
2599                                    // sometimes get a seg fault because
2600                                    // mm is corrupted; strange!?!?!?!?
2601    number nx = number(subst(mm[1], var(1), K[2], var(2), K[3], var(3), K[4]));
2602    number ny = number(subst(mm[2], var(1), K[2], var(2), K[3], var(3), K[4]));
2603    number nz = number(subst(mm[3], var(1), K[2], var(2), K[3], var(3), K[4]));
2604    /* the point (nx, ny, nz) is already a solution;
2605       the following lines will just remove denominators and reduce
2606       numerators in order to return a nice tuple from Z^3 */
2607    bigint nxd = bigint(denominator(absValue(nx)));
2608    bigint nyd = bigint(denominator(absValue(ny)));
2609    bigint nzd = bigint(denominator(absValue(nz)));
2610    nn = nxd * nyd / gcd(nxd, nyd);
2611    nn =  nn * nzd / gcd(nn, nzd);
2612    x = bigint(nx*nn); y = bigint(ny*nn); z = bigint(nz*nn);
2613    nn = gcd(gcd(absValue(x), absValue(y)), absValue(z));
2614    matrix point[1][3];
2615    point[1, 1] = x/nn; point[1, 2] = y/nn; point[1, 3] = z/nn;
2616    export point;
2617    return (R);
2618  }
2619}
2620
2621example
2622{ "EXAMPLE:"; echo=2;
2623ring R = 0, (x,y,z), dp;
2624system("random", 4711);
2625poly p = x^2 + 2*y^2 + 5*z^2 - 4*x*y + 3*x*z + 17*y*z;
2626def S = rationalPointConic(p); // quadratic field extension,
2627                               // minpoly = a^2 - 2
2628testPointConic(p, S);
2629setring R;
2630p = x^2 - 1857669520 * y^2 + 86709575222179747132487270400 * z^2;
2631S = rationalPointConic(p); // same as current basering,
2632                           // no extension needed
2633testPointConic(p, S);
2634}
2635///////////////////////////////////////////////////////////////////////////////
2636proc testParametrization(poly f, def rTT)
2637"USAGE:  testParametrization(f, rTT); f poly, rTT ring
2638ASSUME:  The assumptions on the basering and the polynomial f are as required
2639         by @ref{paraPlaneCurve}. The ring rTT has two variables and contains
2640         an ideal PARA (such as the ring obtained by applying
2641         @ref{paraPlaneCurve} to f).
2642RETURN: int which is 1 if PARA defines a parametrization of the curve
2643        {f=0} and 0, otherwise.
2644THEORY: We compute the polynomial defining the image of PARA
2645        and compare it with f.
2646KEYWORDS: Parametrization, image.
2647EXAMPLE: example testParametrization; shows an example
2648"
2649{
2650  def Roriginal = basering;
2651  setring rTT;
2652  /* begin workaround elimination*/
2653  int k;
2654  list rl = ringlist(rTT);
2655  rl[2] = list("s","t","x","y","z");
2656  rl[3]= list(list("dp",1:5),list("C",0));
2657  def Relim = ring(rl);
2658  setring Relim;
2659  ideal PARA = fetch(rTT,PARA);
2660  ideal JJ;
2661  for(k=1;k<=3;k++)
2662     {
2663       JJ=JJ,var(k+2)-PARA[k];
2664     }
2665  ideal SJJ = std(JJ);
2666  intvec HJJ = hilb(SJJ,1);
2667  ideal J = eliminate(JJ,var(1)*var(2),HJJ);
2668  setring rTT;
2669  /*end workaround elimination*/
2670  rl[2] = list("x","y","z");
2671  rl[3] = list(list("dp",1:3),list("C",0));
2672  def RP2 = ring(rl);
2673  setring RP2;
2674  ideal f = fetch(Roriginal,f);
2675  ideal ftest = imap(Relim,J);
2676  poly g = reduce(f[1],std(ftest));
2677  if(g!=0){return(0)}
2678  g = reduce(ftest[1],std(ideal(f)));
2679  if(g!=0){return(0)}
2680  return (1);
2681}
2682
2683example
2684{ "EXAMPLE:"; echo=2;
2685  ring R = 0,(x,y,z),dp;
2686  poly f = y^8-x^3*(z+x)^5;
2687  def RP1 = paraPlaneCurve(f);
2688  testParametrization(f, RP1);
2689}
2690
2691///////////////////////////////////////////////////////////////////////////////
2692proc testPointConic(poly p, def r)
2693"USAGE:  testPointConic(p, r); p poly, r ring
2694ASSUME:  assumes that p is a homogeneous quadratic polynomial in the
2695        first three ring variables of the current basering;
2696        Assumes that there is a (1x3) matrix named 'point' in r with
2697        entries from the ground field of r.
2698RETURN:  returns 1 iff the point named 'point', residing in r, lies on
2699        the conic given by p; 0 otherwise
2700NOTE:    This method temporarily changes the basering to r. Afterwards,
2701        the basering will be the same as before.
2702EXAMPLE: example testPointConic; shows an example
2703"
2704{
2705 def rOrig = basering;
2706 "conic:", p;
2707 setring r;
2708 string s = "point: " + string(point[1,1]) + ", " + string(point[1,2]);
2709 s = s + ", " + string(point[1,3]);
2710 s;
2711 if (minpoly != 0) { "minpoly:", minpoly; }
2712 poly f = fetch(rOrig, p);
2713 poly g = subst(f, var(1), point[1,1],
2714                   var(2), point[1,2],
2715                   var(3), point[1,3]);
2716 int result = 0; if (g == 0) { result = 1; }
2717 setring rOrig;
2718 return (result);
2719}
2720
2721example
2722{ "EXAMPLE:"; echo=2;
2723 ring R = 0, (x,y,z), dp;
2724 system("random", 4711);
2725 poly p = x^2 + 2*y^2 + 5*z^2 - 4*x*y + 3*x*z + 17*y*z;
2726 def S = rationalPointConic(p);
2727 if (testPointConic(p, S) == 1)
2728 { "point lies on conic"; }
2729 else
2730 { "point does not lie on conic"; }
2731}
2732
2733/////////////////////////////////////////////////////////////////////////////
2734/////////////////////////////////////////////////////////////////////////////
2735/////////////////////////////////////////////////////////////////////////////
2736/*
2737/////////////////////////////////////////////////////////////////////////////
2738/// Further examples for testing the main procedures
2739/// Timings on wawa Sept 29
2740/////////////////////////////////////////////////////////////////////////////
2741LIB"paraplanecurves.lib";
2742// -------------------------------------------------------
2743// Example 1
2744// -------------------------------------------------------
2745ring RR = 0, (x,y,z), dp;
2746poly f = 7*z2+11*y2+13*z*y+17*x2+19*x*y; // conic
2747def RP1 = paraConic(f);
2748setring RP1; PARACONIC;
2749setring RR;
2750RP1 = paraPlaneCurve(f);
2751testParametrization(f,RP1);
2752setring RP1; PARA;
2753kill RR;kill RP1;
2754// -------------------------------------------------------
2755// Example 2
2756// -------------------------------------------------------
2757ring RR = 0, (x,y,z), dp;
2758poly f = y3-x2z;  // cusp at origin
2759adjointIdeal(f,1);
2760adjointIdeal(f,2);
2761def RP1 = paraPlaneCurve(f);  // time 0
2762testParametrization(f,RP1);
2763setring RP1; PARA;
2764kill RR;kill RP1;
2765// -------------------------------------------------------
2766// Example 3
2767// -------------------------------------------------------
2768ring RR = 0, (x,y,z), dp;
2769poly f=(xz-y^2)^2-x*y^3; // 1 sing at origin, 1 cusp, no OMPs
2770adjointIdeal(f,1);
2771adjointIdeal(f,2);
2772def RP1 = paraPlaneCurve(f); // time 0
2773testParametrization(f,RP1);
2774setring RP1; PARA;
2775kill RR;kill RP1;
2776// -------------------------------------------------------
2777// Example 4
2778// -------------------------------------------------------
2779ring RR = 0, (x,y,z), dp;
2780poly f = y5-y4x+4y2x2z-x4z;  // 1 sing at origin, no OMPs, no cusps
2781adjointIdeal(f,1);
2782adjointIdeal(f,2);
2783def RP1 = paraPlaneCurve(f);  // time 0
2784testParametrization(f,RP1);
2785setring RP1; PARA;
2786kill RR;kill RP1;
2787// -------------------------------------------------------
2788// Example 5
2789// -------------------------------------------------------
2790ring RR = 0, (x,y,z), dp;
2791poly f = 259x5-31913x4y+939551x3y2+2871542x2y3+2845801xy4;
2792f = f+914489y5+32068x4z-1884547x3yz-8472623x2y2z-11118524xy3z;
2793f = f-4589347y4z+944585x3z2+8563304x2yz2+16549772xy2z2+9033035y3z2;
2794f = f-2962425x2z3-11214315xyz3-8951744y2z3+2937420xz4+4547571yz4-953955z5;
2795// 6 nodes
2796adjointIdeal(f,1);
2797adjointIdeal(f,2);
2798def RP1 = paraPlaneCurve(f);  // time 0
2799testParametrization(f,RP1);
2800setring RP1; PARA;
2801kill RR;kill RP1;
2802// -------------------------------------------------------
2803// Example 7
2804// -------------------------------------------------------
2805ring RR = 0, (x,y,z), dp;
2806poly f = y^8-x^3*(z+x)^5;  // 1 sing at origin, 1 further sing, no OMPs,
2807                           // no cusps
2808adjointIdeal(f,1);
2809adjointIdeal(f,2);
2810def RP1 = paraPlaneCurve(f);  // time 0
2811testParametrization(f,RP1);
2812setring RP1; PARA;
2813kill RR;kill RP1;
2814// -------------------------------------------------------
2815// Example 8
2816// -------------------------------------------------------
2817ring RR = 0, (x,y,z), dp;
2818poly f = 11y7+7y6x+8y5x2-3y4x3-10y3x4-10y2x5-x7-33y6-29y5x-13y4x2+26y3x3;
2819f = f+30y2x4+10yx5+3x6+33y5+37y4x-8y3x2-33y2x3-20yx4-3x5-11y4-15y3x;
2820f = f+13y2x2+10yx3+x4; // 3 OMPs of mult 3, 1 OMP of mult 4
2821f = homog(f,z);
2822adjointIdeal(f,1);
2823adjointIdeal(f,2);
2824def RP1 = paraPlaneCurve(f);  // time 0
2825testParametrization(f,RP1);
2826setring RP1; PARA;
2827kill RR;kill RP1;
2828// -------------------------------------------------------
2829// Example 9
2830// -------------------------------------------------------
2831ring RR = 0, (x,y,z), dp;
2832poly f = y^8-x^3*(z+x)^5;  // 1 sing at origin, 1 further sing, no OMPs,
2833                           // no cusps
2834adjointIdeal(f,1);
2835adjointIdeal(f,2);
2836def RP1 = paraPlaneCurve(f);  // time 0
2837testParametrization(f,RP1);
2838setring RP1; PARA;
2839kill RR;kill RP1;
2840// -------------------------------------------------------
2841// Example 10
2842// -------------------------------------------------------
2843ring SS = 0, (u,v,z), dp;
2844poly f = u^4-14*u^2*v^2+v^4+8*u^2*v*z+8*v^3*z; // 1 OMP of mult 3 at orgin
2845adjointIdeal(f,1);
2846adjointIdeal(f,2);
2847def RP1 = paraPlaneCurve(f);  // time 0
2848testParametrization(f,RP1);
2849setring RP1; PARA;
2850kill SS;kill RP1;
2851// -------------------------------------------------------
2852// Example 11
2853// -------------------------------------------------------
2854ring SS = 0, (u,v,z), dp;
2855poly f = 14440*u^5-16227*u^4*v+10812*u^3*v^2-13533*u^2*v^3+3610*u*v^4;
2856f = f+1805*v^5+14440*u^4*z-18032*u^3*v*z+16218*u^2*v^2*z-12626*u*v^3*z;
2857f = f+3610*v^4*z+3610*u^3*z^2-4508*u^2*v*z^2+5406*u*v^2*z^2-2703*v^3*z^2;
2858// 1 OMP of mult 3 at origin, 2 nodes
2859adjointIdeal(f,1);
2860adjointIdeal(f,2);
2861def RP1 = paraPlaneCurve(f);  // time 0
2862testParametrization(f,RP1);
2863setring RP1; PARA;
2864kill SS;kill RP1;
2865// -------------------------------------------------------
2866// Example 12
2867// -------------------------------------------------------
2868ring SS = 0, (u,v,z), dp;
2869poly f = u^6+3*u^4*v^2+3*u^2*v^4+v^6-4*u^4*z^2-34*u^3*v*z^2-7*u^2*v^2*z^2;
2870f = f+12*u*v^3*z^2+6*v^4*z^2+36*u^2*z^4+36*u*v*z^4+9*v^2*z^4;
2871// needs field extension *** 6 nodes, 2 cusps, 1 sing at 0
2872adjointIdeal(f,1);
2873adjointIdeal(f,2);
2874def RP1 = paraPlaneCurve(f);  // time 0
2875testParametrization(f,RP1);
2876setring RP1; PARA;
2877kill SS;kill RP1;
2878// -------------------------------------------------------
2879// Example 13
2880// -------------------------------------------------------
2881ring SS = 0, (u,v,z), dp;
2882poly f = -24135/322*u^6-532037/6440*u^5*v+139459/560*u^4*v^2;
2883f = f-1464887/12880*u^3*v^3+72187/25760*u^2*v^4+9/8*u*v^5+1/8*v^6;
2884f = f-403511/3220*u^5*z-40817/920*u^4*v*z+10059/80*u^3*v^2*z;
2885f = f-35445/1288*u^2*v^3*z+19/4*u*v^4*z+3/4*v^5*z-20743/805*u^4*z^2;
2886f = f+126379/3220*u^3*v*z^2-423417/6440*u^2*v^2*z^2+11/2*u*v^3*z^2;
2887f = f+3/2*v^4*z^2+3443/140*u^3*z^3+u^2*v*z^3+u*v^2*z^3+v^3*z^3;
2888// 2 OMPs of mult 3 (1 at origin), 4 nodes
2889adjointIdeal(f,1);
2890adjointIdeal(f,2);
2891def RP1 = paraPlaneCurve(f);  // time 14
2892testParametrization(f,RP1);
2893setring RP1; PARA;
2894kill SS;kill RP1;
2895// -------------------------------------------------------
2896// Example 14
2897// -------------------------------------------------------
2898ring SS = 0, (u,v,z), dp;
2899poly f =
29002*u^7+u^6*v+3*u^5*v^2+u^4*v^3+2*u^3*v^4+u^2*v^5+2*u*v^6+v^7
2901-7780247/995328*u^6*z-78641/9216*u^5*v*z-10892131/995328*u^4*v^2*z
2902-329821/31104*u^3*v^3*z-953807/331776*u^2*v^4*z-712429/248832*u*v^5*z
2903+1537741/331776*v^6*z+2340431/248832*u^5*z^2+5154337/248832*u^4*v*z^2
2904+658981/41472*u^3*v^2*z^2+1737757/124416*u^2*v^3*z^2
2905-1234733/248832*u*v^4*z^2-1328329/82944*v^5*z^2-818747/248832*u^4*z^3
2906-1822879/124416*u^3*v*z^3-415337/31104*u^2*v^2*z^3
2907+1002655/124416*u*v^3*z^3+849025/82944*v^4*z^3;
2908// 3 OMPs of mult 3, 1 OMP of mult 4 at origin
2909adjointIdeal(f,2);
2910def RP1 = paraPlaneCurve(f);  // time 1
2911testParametrization(f,RP1);
2912setring RP1; PARA;
2913kill SS;kill RP1;
2914// -------------------------------------------------------
2915// Example 15
2916// -------------------------------------------------------
2917ring SS = 0, (u,v,z), dp;
2918poly f = 590819418867856650536224u7-147693905508217596067968u6v;
2919f = f+229117518934972047619978u5v2-174050799674982973889542u4v3;
2920f = f-92645796479789150855110u3v4-65477418713685583062704u2v5;
2921f = f+4529961835917468460168uv6+7715404057796585983136v7;
2922f = f-413640780091141905428104u6z+571836835577486968144618u5vz;
2923f = f-551807810327826605739444u4v2z-488556410340789283359926u3v3z;
2924f = f-473466023008413178155962u2v4z+48556741573432247323608uv5z;
2925f = f+77647371229172269259528v6z+340450118906560552282893u5z2;
2926f = f-433598825064368371610344u4vz2-937281070591684636591672u3v2z2;
2927f = f-1388949843915129934647751u2v3z2+204081793110898617103998uv4z2;
2928f = f+335789953068251652554308v5z2+6485661002496681852577u4z3;
2929f = f-772700266516318390630202u3vz3-2068348417248100329533330u2v2z3;
2930f = f+440320154612359641806108uv3z3+808932515589210854581618v4z3;
2931f = f-229384307132237615286548u3z4-1564303565658228216055227u2vz4;
2932f = f+520778334468674798322974uv2z4+1172483905704993294097655v3z4;
2933f = f-480789741398016816562100u2z5+322662751598958620410786uvz5;
2934f = f+1022525576391791616258310v2z5+82293493608853837667471uz6;
2935f = f+496839109904761426785889vz6+103766136235628614937587z7; // 15 nodes
2936adjointIdeal(f,2);
2937def RP1 = paraPlaneCurve(f);  // time 72
2938testParametrization(f,RP1);
2939setring RP1; PARA;
2940kill SS;kill RP1;
2941
2942// -------------------------------------------------------
2943// Example 16
2944// -------------------------------------------------------
2945ring SS = 0, (u,v,z), dp;
2946poly f = 25*u^8+184*u^7*v+518*u^6*v^2+720*u^5*v^3+576*u^4*v^4+282*u^3*v^5;
2947f = f+84*u^2*v^6+14*u*v^7+v^8+244*u^7*z+1326*u^6*v*z+2646*u^5*v^2*z;
2948f = f+2706*u^4*v^3*z+1590*u^3*v^4*z+546*u^2*v^5*z+102*u*v^6*z+8*v^7*z;
2949f = f+854*u^6*z^2+3252*u^5*v*z^2+4770*u^4*v^2*z^2+3582*u^3*v^3*z^2;
2950f = f+1476*u^2*v^4*z^2+318*u*v^5*z^2+28*v^6*z^2+1338*u^5*z^3+3740*u^4*v*z^3;
2951f = f+4030*u^3*v^2*z^3+2124*u^2*v^3*z^3+550*u*v^4*z^3+56*v^5*z^3+1101*u^4*z^4;
2952f = f+2264*u^3*v*z^4+1716*u^2*v^2*z^4+570*u*v^3*z^4+70*v^4*z^4+508*u^3*z^5;
2953f = f+738*u^2*v*z^5+354*u*v^2*z^5+56*v^3*z^5+132*u^2*z^6+122*u*v*z^6;
2954f = f+28*v^2*z^6+18*u*z^7+8*v*z^7+z^8; // 3 nodes, 1 sing
2955adjointIdeal(f,1);
2956adjointIdeal(f,2);
2957def RP1 = paraPlaneCurve(f);  // time 20
2958testParametrization(f,RP1);
2959setring RP1; PARA;
2960kill SS;kill RP1;
2961// -------------------------------------------------------
2962// Example 17
2963// -------------------------------------------------------
2964ring SS = 0, (u,v,z), dp;
2965poly f = -2*u*v^4*z^4+u^4*v^5+12*u^4*v^3*z^2+12*u^2*v^4*z^3-u^3*v*z^5;
2966f = f+11*u^3*v^2*z^4-21*u^3*v^3*z^3-4*u^4*v*z^4+2*u^4*v^2*z^3-6*u^4*v^4*z;
2967f = f+u^5*z^4-3*u^5*v^2*z^2+u^5*v^3*z-3*u*v^5*z^3-2*u^2*v^3*z^4+u^3*v^4*z^2;
2968f = f+v^5*z^4; // 2 OMPs of mult 4, 1 OMP of mult 5, 1 sing at origin
2969f = subst(f,z,u+v+z);
2970adjointIdeal(f,2);
2971def RP1 = paraPlaneCurve(f);  // time 5
2972testParametrization(f,RP1);
2973setring RP1; PARA;
2974kill SS;kill RP1;
2975// -------------------------------------------------------
2976// Example 18
2977// -------------------------------------------------------
2978ring SS = 0, (u,v,z), dp;
2979poly f = u^5*v^5+21*u^5*v^4*z-36*u^4*v^5*z-19*u^5*v^3*z^2+12*u^4*v^4*z^2;
2980f = f+57*u^3*v^5*z^2+u^5*v^2*z^3+u^4*v^3*z^3-53*u^3*v^4*z^3-19*u^2*v^5*z^3;
2981f = f+u^5*v*z^4+43*u^3*v^3*z^4+u*v^5*z^4+u^5*z^5-15*u^3*v^2*z^5+u^2*v^3*z^5;
2982f = f+u*v^4*z^5+v^5*z^5; // 1 OMP of mult 4, 3 OMPs of mult 5 (1 at origin)
2983adjointIdeal(f,2);
2984def RP1 = paraPlaneCurve(f);  // time 8
2985testParametrization(f,RP1);
2986setring RP1; PARA;
2987kill SS;kill RP1;
2988// -------------------------------------------------------
2989// Example 19
2990// -------------------------------------------------------
2991ring SS = 0, (u,v,z), dp;
2992poly f = u^10+6*u^9*v-30*u^7*v^3-15*u^6*v^4+u^5*v^5+u^4*v^6+6*u^3*v^7;
2993f = f+u^2*v^8+7*u*v^9+v^10+5*u^9*z+24*u^8*v*z-30*u^7*v^2*z-120*u^6*v^3*z;
2994f = f-43*u^5*v^4*z+5*u^4*v^5*z+20*u^3*v^6*z+10*u^2*v^7*z+29*u*v^8*z+5*v^9*z;
2995f = f+10*u^8*z^2+36*u^7*v*z^2-105*u^6*v^2*z^2-179*u^5*v^3*z^2-38*u^4*v^4*z^2;
2996f = f+25*u^3*v^5*z^2+25*u^2*v^6*z^2+46*u*v^7*z^2+10*v^8*z^2+10*u^7*z^3;
2997f = f+24*u^6*v*z^3-135*u^5*v^2*z^3-117*u^4*v^3*z^3-u^3*v^4*z^3+25*u^2*v^5*z^3;
2998f = f+34*u*v^6*z^3+10*v^7*z^3+5*u^6*z^4+6*u^5*v*z^4-75*u^4*v^2*z^4;
2999f = f-27*u^3*v^3*z^4+10*u^2*v^4*z^4+11*u*v^5*z^4+5*v^6*z^4+u^5*z^5;
3000f = f-15*u^3*v^2*z^5+u^2*v^3*z^5+u*v^4*z^5+v^5*z^5;
3001// 1 OMP of mult 4, 3 OMPs of mult 5 (1 at origin)
3002adjointIdeal(f,2);
3003def RP1 = paraPlaneCurve(f);  // time 2 // see Ex. 18
3004testParametrization(f,RP1);
3005setring RP1; PARA;
3006kill SS;kill RP1;
3007// -------------------------------------------------------
3008// Example 20
3009// -------------------------------------------------------
3010ring R = 0, (x,y,z), dp;
3011system("random", 4711);
3012poly p = x^2 + 2*y^2 + 5*z^2 - 4*x*y + 3*x*z + 17*y*z;
3013def S = rationalPointConic(p); // quadratic field extension,
3014                              // minpoly = a^2 - 2
3015if (testPointConic(p, S) == 1)
3016{ "point lies on conic"; }
3017else
3018{ "point does not lie on conic"; }
3019kill R;kill S;
3020// -------------------------------------------------------
3021// Example 21
3022// -------------------------------------------------------
3023ring R = 0, (x,y,z), dp;
3024system("random", 4711);
3025poly p = x^2 - 1857669520 * y^2 + 86709575222179747132487270400 * z^2;
3026def S = rationalPointConic(p); // same as current basering,
3027                              // no extension needed
3028if (testPointConic(p, S) == 1)
3029{ "point lies on conic"; }
3030else
3031{ "point does not lie on conic"; }
3032kill R;kill S;
3033// -------------------------------------------------------
3034// Example 21
3035// -------------------------------------------------------
3036ring RR = 0, (x,y,z), dp;
3037poly f = -1965466244509920x5y+34871245546721380061760x4y2;
3038f = f+104613747941595046117320x3y3+113331564241941002407560x2y4;
3039f = f+52306876673313609259800xy5+8717812860780028397880y6;
3040f = f+1040297748510024x5z+4468147845634872x4yz;
3041f = f-22398508728211453743258x3y2z-33223996581074443306854x2y3z;
3042f = f-10638598235041298082366xy4z+186886189971594356382y5z;
3043f = f-1385078844909312x4z2-34893092731637052532683x3yz2;
3044f = f-98591463214095439056609x2y2z2-92339459334829609336485xy3z2;
3045f = f-24923289542522905755711y4z2+472440640471377x3z3;
3046f = f+33821511925664516716011x2yz3+49745237303968344397437xy2z3;
3047f = f+11040465960074786720475y3z3+8728735735878837099404x2z4;
3048f = f+17676785754519678518537xyz4+17935885079051421934609y2z4;
3049f = f-11314701999743172607075xz5-16164284825803158969425yz5;
3050f = f+3666695988537425618750z6;
3051// 4 nodes, 1 OMP of mult 4
3052adjointIdeal(f,2);
3053kill RR;
3054*/
Note: See TracBrowser for help on using the repository browser.