source: git/Singular/LIB/paraplanecurves.lib @ 6b5e56

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