source: git/Singular/LIB/paraplanecurves.lib @ bee06d

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