source: git/Singular/LIB/paraplanecurves.lib @ 1a96a4

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