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

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